ATMOS_PHY_BL_MYNN_tendency calculate tendency by the vertical eddy viscosity.
725 hydrometeor_lhv => atmos_hydrometeor_lhv, &
728 matrix_solver_tridiagonal
730 file_history_query, &
734 integer,
intent(in) :: KA, KS, KE
735 integer,
intent(in) :: IA, IS, IE
736 integer,
intent(in) :: JA, JS, JE
738 real(RP),
intent(in) :: DENS (KA,IA,JA)
739 real(RP),
intent(in) :: U (KA,IA,JA)
740 real(RP),
intent(in) :: V (KA,IA,JA)
741 real(RP),
intent(in) :: W (KA,IA,JA)
742 real(RP),
intent(in) :: POTT (KA,IA,JA)
743 real(RP),
intent(in) :: PROG (KA,IA,JA,ATMOS_PHY_BL_MYNN_ntracer)
744 real(RP),
intent(in) :: PRES (KA,IA,JA)
745 real(RP),
intent(in) :: EXNER (KA,IA,JA)
746 real(RP),
intent(in) :: N2 (KA,IA,JA)
747 real(RP),
intent(in) :: QDRY (KA,IA,JA)
748 real(RP),
intent(in) :: QV (KA,IA,JA)
749 real(RP),
intent(in) :: Qw (KA,IA,JA)
750 real(RP),
intent(in) :: POTL (KA,IA,JA)
751 real(RP),
intent(in) :: POTV (KA,IA,JA)
752 real(RP),
intent(in) :: SFC_DENS( IA,JA)
753 real(RP),
intent(in) :: SFLX_MU ( IA,JA)
754 real(RP),
intent(in) :: SFLX_MV ( IA,JA)
755 real(RP),
intent(in) :: SFLX_SH ( IA,JA)
756 real(RP),
intent(in) :: SFLX_QV ( IA,JA)
757 real(RP),
intent(in) :: us ( IA,JA)
758 real(RP),
intent(in) :: ts ( IA,JA)
759 real(RP),
intent(in) :: qs ( IA,JA)
760 real(RP),
intent(in) :: RLmo ( IA,JA)
762 real(RP),
intent(in) :: frac_land(IA,JA)
764 real(RP),
intent(in) :: CZ( KA,IA,JA)
765 real(RP),
intent(in) :: FZ(0:KA,IA,JA)
766 real(RP),
intent(in) :: F2H(KA,2,IA,JA)
767 real(DP),
intent(in) :: dt_DP
769 character(len=*),
intent(in) :: BULKFLUX_type
771 real(RP),
intent(out) :: RHOU_t (KA,IA,JA)
772 real(RP),
intent(out) :: RHOV_t (KA,IA,JA)
773 real(RP),
intent(out) :: RHOT_t (KA,IA,JA)
774 real(RP),
intent(out) :: RHOQV_t(KA,IA,JA)
775 real(RP),
intent(out) :: RPROG_t(KA,IA,JA,ATMOS_PHY_BL_MYNN_ntracer)
776 real(RP),
intent(out) :: Nu (KA,IA,JA)
777 real(RP),
intent(out) :: Kh (KA,IA,JA)
778 real(RP),
intent(out) :: Qlp (KA,IA,JA)
779 real(RP),
intent(out) :: cldfrac(KA,IA,JA)
780 real(RP),
intent(out) :: Zi (IA,JA)
781 real(RP),
intent(out) :: SFLX_BUOY (IA,JA)
786 real(RP) :: Ri (KA,IA,JA)
787 real(RP) :: Pr (KA,IA,JA)
788 real(RP) :: prod (KA,IA,JA)
789 real(RP) :: diss (KA,IA,JA)
790 real(RP) :: dudz2(KA,IA,JA)
791 real(RP) :: l (KA,IA,JA)
792 real(RP) :: rho_h(KA)
794 real(RP) :: flxU(0:KA,IA,JA)
795 real(RP) :: flxV(0:KA,IA,JA)
796 real(RP) :: flxT(0:KA,IA,JA)
797 real(RP) :: flxQ(0:KA,IA,JA)
799 real(RP) :: flxU2(0:KA,IA,JA)
800 real(RP) :: flxV2(0:KA,IA,JA)
801 real(RP) :: flxT2(0:KA,IA,JA)
802 real(RP) :: flxQ2(0:KA,IA,JA)
805 real(RP) :: RHONu (KA)
806 real(RP) :: RHOKh (KA)
807 real(RP) :: N2_new(KA)
809 real(RP) :: sm25 (KA)
810 real(RP) :: sh25 (KA)
820 real(RP) :: prod_t(KA)
821 real(RP) :: prod_q(KA)
822 real(RP) :: prod_c(KA)
823 real(RP) :: diss_p(KA)
824 real(RP) :: dtldz(KA)
825 real(RP) :: dqwdz(KA)
827 real(RP) :: shpgh(KA)
828 real(RP) :: gammat (KA)
829 real(RP) :: gammaq (KA)
830 real(RP) :: rlqsm_h(KA)
832 real(RP) :: flx(0:KA)
838 real(RP) :: phi_N(KA)
839 real(RP) :: dummy(KA)
849 logical :: mynn_level3
857 real(RP) :: tflux(0:KA)
858 real(RP) :: qflux(0:KA)
859 real(RP) :: uflux(0:KA)
860 real(RP) :: vflux(0:KA)
861 real(RP) :: eflux(0:KA)
872 real(RP) :: PTLV (KA)
873 real(RP) :: Nu_f (KA)
874 real(RP) :: Kh_f (KA)
875 real(RP) :: q2_2 (KA)
877 real(RP) :: betat(KA)
878 real(RP) :: betaq(KA)
880 real(RP) :: f_smp (KA)
881 real(RP) :: f_shpgh(KA)
882 real(RP) :: f_gamma(KA)
883 real(RP) :: tltv25 (KA)
884 real(RP) :: qwtv25 (KA)
885 real(RP) :: tvsq25 (KA)
886 real(RP) :: tvsq_up(KA)
887 real(RP) :: tvsq_lo(KA)
888 real(RP) :: dtsq (KA)
889 real(RP) :: dqsq (KA)
890 real(RP) :: dcov (KA)
892 real(RP) :: mflux(0:KA)
894 real(RP) :: Uh(KA), Vh(KA), Wh(KA)
895 real(RP) :: qw2(KA), qh(KA)
896 real(RP) :: tlh(KA), tvh(KA)
897 real(RP) :: eh(KA), dh(KA)
898 real(RP) :: TEML(KA), LHVL(KA), psat(KA)
901 real(RP) :: work(KA,4)
908 log_progress(*)
"atmosphere / physics / pbl / MYNN"
910 mynn_level3 = ( atmos_phy_bl_mynn_level ==
"3" )
912 select case ( bulkflux_type )
920 log_error(
"ATMOS_PHY_BL_MYNN_tendency",*)
"BULKFLUX_type is invalid: ", trim(bulkflux_type)
994 z(k) = cz(k,i,j) - fz(ks-1,i,j)
1000 if ( atmos_phy_bl_mynn_pbl_max >= z(k) )
then
1008 cdz(k) = fz(k,i,j) - fz(k-1,i,j)
1011 fdz(k) = cz(k+1,i,j) - cz(k,i,j)
1015 tke(k) = max( prog(k,i,j,
i_tke), atmos_phy_bl_mynn_tke_min )
1018 if ( mynn_level3 )
then
1020 tsq(k) = max( prog(k,i,j,i_tsq), 0.0_rp )
1021 qsq(k) = max( prog(k,i,j,i_qsq), 0.0_rp )
1022 cov(k) = prog(k,i,j,i_cov)
1023 cov(k) = sign( min( abs(cov(k)), sqrt(tsq(k) * qsq(k))), cov(k) )
1027 cptot = cpdry + sflx_qv(i,j) * ( cp_vapor - cpdry )
1028 sflx_pt = sflx_sh(i,j) / ( cptot * exner(ks,i,j) )
1030 call mynn_main( ka, ks, ke_pbl, &
1032 tke(:), tsq(:), qsq(:), cov(:), &
1033 q(:), l(:,i,j), lq(:), &
1034 nu(:,i,j), rhonu(:), kh(:,i,j), rhokh(:), pr(:,i,j), &
1035 prod(:,i,j), prod_t(:), prod_q(:), prod_c(:), &
1036 diss(:,i,j), diss_p(:), &
1037 sm25(:), smp(:), sh25(:), shpgh(:), &
1038 gammat(:), gammaq(:), &
1039 dudz2(:,i,j), n2_new(:), ri(:,i,j), &
1040 dtldz(:), dqwdz(:), &
1042 uflux(:), vflux(:), tflux(:), qflux(:), eflux(:), &
1043 qlp(:,i,j), cldfrac(:,i,j), zi(i,j), sflx_buoy(i,j), &
1044 u(:,i,j), v(:,i,j), w(:,i,j), &
1045 dens(:,i,j), pres(:,i,j), &
1046 pott(:,i,j), potl(:,i,j), potv(:,i,j), &
1047 qw(:,i,j), n2(:,i,j), &
1048 exner(:,i,j), qdry(:,i,j), &
1049 sflx_pt, sflx_sh(i,j), sflx_qv(i,j), sfc_dens(i,j), &
1050 rlmo(i,j), us(i,j), ts(i,j), qs(i,j), &
1051 z(:), cdz(:), fdz(:), f2h(:,:,i,j), &
1055 mynn_level3, .false., &
1059 ptlv(:), nu_f(:), kh_f(:), q2_2(:), ac(:), &
1060 betat(:), betaq(:), &
1061 flx(:), a(:), b(:), c(:), d(:), &
1062 f_smp(:), f_shpgh(:), f_gamma(:), &
1063 tltv25(:), qwtv25(:), tvsq25(:), &
1064 tvsq_up(:), tvsq_lo(:), dtsq(:), dqsq(:), dcov(:), &
1066 uh(:), vh(:), wh(:), qw2(:), qh(:), &
1067 tlh(:), tvh(:), eh(:), dh(:), &
1068 teml(:), lhvl(:), psat(:) )
1074 flx(ke_pbl) = 0.0_rp
1077 rlqsm_h(k) = rho_h(k) * ( f2h(k,1,i,j) * lq(k+1) * smp(k+1) + f2h(k,2,i,j) * lq(k) * smp(k) )
1082 if ( atmos_phy_bl_mynn_mf )
then
1089 flx(k) = - rlqsm_h(k) * ( u(k+1,i,j) - u(k,i,j) ) / fdz(k)
1093 sf_t = sflx_mu(i,j) / cdz(ks)
1094 d(ks) = ( u(ks,i,j) * dens(ks,i,j) + dt * ( sf_t - flx(ks) / cdz(ks) ) ) / rho(ks)
1096 d(k) = u(k,i,j) - dt * ( flx(k) - flx(k-1) ) / ( cdz(k) * rho(k) )
1100 ap = - dt * rhonu(k) / fdz(k)
1101 a(k) = ap / ( rho(k) * cdz(k) )
1104 b(k) = - a(k) + 1.0_rp
1106 b(k) = - a(k) + dt * rhonu(k-1) / ( fdz(k-1) * rho(k) * cdz(k) ) + 1.0_rp
1109 b(k) = - a(k) - c(k) + 1.0_rp
1111 c(k+1) = ap / ( rho(k+1) * cdz(k+1) )
1114 b(ke_pbl) = - c(ke_pbl) + 1.0_rp
1116 call matrix_solver_tridiagonal( &
1121 a(:), b(:), c(:), d(:), &
1125 phi_n(ks:ke_pbl) = dummy(ks:ke_pbl)
1126 rhou_t(ks,i,j) = ( phi_n(ks) * rho(ks) - u(ks,i,j) * dens(ks,i,j) ) / dt - sf_t
1128 rhou_t(k,i,j) = ( phi_n(k) - u(k,i,j) ) * rho(k) / dt
1131 rhou_t(k,i,j) = 0.0_rp
1133 flxu(ks-1,i,j) = 0.0_rp
1134 flxu2(ks-1,i,j) = 0.0_rp
1136 flxu(k,i,j) = flx(k) &
1137 - rhonu(k) * ( phi_n(k+1) - phi_n(k) ) / fdz(k)
1138 flxu2(k,i,j) = flx(k)
1141 flxu(k,i,j) = 0.0_rp
1142 flxu2(k,i,j) = 0.0_rp
1148 if ( atmos_phy_bl_mynn_mf )
then
1155 flx(k) = - rlqsm_h(k) * ( v(k+1,i,j) - v(k,i,j) ) / fdz(k)
1159 sf_t = sflx_mv(i,j) / cdz(ks)
1160 d(ks) = ( v(ks,i,j) * dens(ks,i,j) + dt * ( sf_t - flx(ks) / cdz(ks) ) ) / rho(ks)
1162 d(k) = v(k,i,j) - dt * ( flx(k) - flx(k-1) ) / ( cdz(k) * rho(k) )
1166 call matrix_solver_tridiagonal( &
1171 a(:), b(:), c(:), d(:), &
1174 rhov_t(ks,i,j) = ( phi_n(ks) * rho(ks) - v(ks,i,j) * dens(ks,i,j) ) / dt - sf_t
1176 rhov_t(k,i,j) = ( phi_n(k) - v(k,i,j) ) * rho(k) / dt
1179 rhov_t(k,i,j) = 0.0_rp
1181 flxv(ks-1,i,j) = 0.0_rp
1182 flxv2(ks-1,i,j) = 0.0_rp
1184 flxv(k,i,j) = flx(k) &
1185 - rhonu(k) * ( phi_n(k+1) - phi_n(k) ) / fdz(k)
1186 flxv2(k,i,j) = flx(k)
1189 flxv(k,i,j) = 0.0_rp
1190 flxv2(k,i,j) = 0.0_rp
1196 if ( atmos_phy_bl_mynn_mf )
then
1203 flx(k) = - ( f2h(k,1,i,j) * lq(k+1) * gammat(k+1) + f2h(k,2,i,j) * lq(k) * gammat(k) ) &
1208 sf_t = sflx_pt / cdz(ks)
1210 d(ks) = pott(ks,i,j) + dt * ( sf_t - flx(ks) / cdz(ks) ) / rho(ks)
1212 d(k) = pott(k,i,j) - dt * ( flx(k) - flx(k-1) ) / ( cdz(k) * rho(k) )
1217 ap = - dt * rhokh(k) / fdz(k)
1218 a(k) = ap / ( rho(k) * cdz(k) )
1221 b(k) = - a(k) + 1.0_rp
1223 b(k) = - a(k) + dt * rhokh(k-1) / ( fdz(k-1) * rho(k) * cdz(k) ) + 1.0_rp
1226 b(k) = - a(k) - c(k) + 1.0_rp
1228 c(k+1) = ap / ( rho(k+1) * cdz(k+1) )
1231 b(ke_pbl) = - c(ke_pbl) + 1.0_rp
1233 call matrix_solver_tridiagonal( &
1238 a(:), b(:), c(:), d(:), &
1241 rhot_t(ks,i,j) = ( phi_n(ks) - pott(ks,i,j) ) * rho(ks) / dt - sf_t
1243 rhot_t(k,i,j) = ( phi_n(k) - pott(k,i,j) ) * rho(k) / dt
1246 rhot_t(k,i,j) = 0.0_rp
1248 flxt(ks-1,i,j) = 0.0_rp
1249 flxt2(ks-1,i,j) = 0.0_rp
1251 flxt(k,i,j) = flx(k) &
1252 - rhokh(k) * ( phi_n(k+1) - phi_n(k) ) / fdz(k)
1253 flxt2(k,i,j) = flx(k)
1256 flxt(k,i,j) = 0.0_rp
1257 flxt2(k,i,j) = 0.0_rp
1261 if ( flxt(ks,i,j) > 1e-4_rp )
then
1262 fmin = flxt(ks,i,j) / rho_h(ks)
1264 do k = ks+1, ke_pbl-2
1265 tmp = ( flxt(k-1,i,j) + flxt(k,i,j) + flxt(k+1,i,j) ) &
1266 / ( rho_h(k-1) + rho_h(k) + rho_h(k+1) )
1267 if ( fmin < 0.0_rp .and. tmp > fmin )
exit
1268 if ( tmp < fmin )
then
1274 zi(i,j) = fz(kmin,i,j) - fz(ks-1,i,j)
1278 if ( atmos_phy_bl_mynn_mf )
then
1285 flx(k) = - ( f2h(k,1,i,j) * lq(k+1) * gammaq(k+1) + f2h(k,2,i,j) * lq(k) * gammaq(k) ) &
1290 sf_t = sflx_qv(i,j) / cdz(ks)
1291 d(ks) = ( qv(ks,i,j) * dens(ks,i,j) + dt * ( sf_t - flx(ks) / cdz(ks) ) ) / rho(ks)
1293 d(k) = qv(k,i,j) - dt * ( flx(k) - flx(k-1) ) / ( cdz(k) * rho(k) )
1296 call matrix_solver_tridiagonal( &
1301 a(:), b(:), c(:), d(:), &
1304 rhoqv_t(ks,i,j) = ( phi_n(ks) * rho(ks) - qv(ks,i,j) * dens(ks,i,j) ) / dt - sf_t
1306 rhoqv_t(k,i,j) = ( phi_n(k) - qv(k,i,j) ) * rho(k) / dt
1309 rhoqv_t(k,i,j) = 0.0_rp
1311 flxq(ks-1,i,j) = 0.0_rp
1312 flxq2(ks-1,i,j) = 0.0_rp
1314 flxq(k,i,j) = flx(k) &
1315 - rhokh(k) * ( phi_n(k+1) - phi_n(k) ) / fdz(k)
1316 flxq2(k,i,j) = flx(k)
1319 flxq(k,i,j) = 0.0_rp
1320 flxq2(k,i,j) = 0.0_rp
1336 d(k) = tke(k) * dens(k,i,j) / rho(k) + dt * ( - ( flx(k) - flx(k-1) ) / ( cdz(k) * rho(k) ) + prod(k,i,j) )
1342 ap = - dt * atmos_phy_bl_mynn_sq_fact * rhonu(k) / fdz(k)
1343 a(k) = ap / ( rho(k) * cdz(k) )
1346 b(k) = - a(k) + 1.0_rp - diss(k,i,j) * dt
1348 b(k) = - a(k) + dt * atmos_phy_bl_mynn_sq_fact * rhonu(k-1) / ( fdz(k-1) * rho(k) * cdz(k) ) + 1.0_rp - diss(k,i,j) * dt
1351 b(k) = - a(k) - c(k) + 1.0_rp - diss(k,i,j) * dt
1353 c(k+1) = ap / ( rho(k+1) * cdz(k+1) )
1357 call matrix_solver_tridiagonal( &
1362 a(:), b(:), c(:), d(:), &
1366 phi_n(k) = max( phi_n(k), atmos_phy_bl_mynn_tke_min )
1370 diss(k,i,j) = diss(k,i,j) * phi_n(k)
1371 rprog_t(k,i,j,
i_tke) = ( phi_n(k) * rho(k) - prog(k,i,j,
i_tke) * dens(k,i,j) ) / dt
1375 rprog_t(k,i,j,
i_tke) = - atmos_phy_bl_mynn_dump_coef * prog(k,i,j,
i_tke) * dens(k,i,j) / dt
1376 diss(k,i,j) = 0.0_rp
1377 prod(k,i,j) = 0.0_rp
1381 if ( mynn_level3 )
then
1386 diss_p(k) = dt * 2.0_rp * q(k) / ( b2 * l(k,i,j) )
1389 d(k) = tsq(k) * dens(k,i,j) / rho(k) + dt * prod_t(k)
1394 ap = - dt * rhonu(k) / fdz(k)
1395 a(k) = ap / ( rho(k) * cdz(k) )
1398 b(k) = - a(k) + diss_p(k)
1400 b(k) = - a(k) + dt * rhonu(k-1) / ( fdz(k-1) * rho(k) * cdz(k) ) + diss_p(k)
1403 b(k) = - a(k) - c(k) + 1.0_rp + diss_p(k)
1405 c(k+1) = ap / ( rho(k+1) * cdz(k+1) )
1408 b(ke_pbl) = - c(ke_pbl) + 1.0_rp + diss_p(ke_pbl)
1410 call matrix_solver_tridiagonal( &
1415 a(:), b(:), c(:), d(:), &
1419 tsq(k) = max( tsq(k), 0.0_rp )
1423 rprog_t(k,i,j,i_tsq) = ( tsq(k) * rho(k) - prog(k,i,j,i_tsq) * dens(k,i,j) ) / dt
1427 rprog_t(k,i,j,i_tsq) = - atmos_phy_bl_mynn_dump_coef * prog(k,i,j,i_tsq) * dens(k,i,j) / dt
1434 d(k) = qsq(k) * dens(k,i,j) / rho(k) + dt * prod_q(k)
1438 call matrix_solver_tridiagonal( &
1443 a(:), b(:), c(:), d(:), &
1447 qsq(k) = max( qsq(k), 0.0_rp )
1451 rprog_t(k,i,j,i_qsq) = ( qsq(k) * rho(k) - prog(k,i,j,i_qsq) * dens(k,i,j) ) / dt
1455 rprog_t(k,i,j,i_qsq) = - atmos_phy_bl_mynn_dump_coef * prog(k,i,j,i_qsq) * dens(k,i,j) / dt
1462 d(k) = cov(k) * dens(k,i,j) / rho(k) + dt * prod_c(k)
1467 call matrix_solver_tridiagonal( &
1472 a(:), b(:), c(:), d(:), &
1476 cov(k) = sign( min( abs(cov(k)), sqrt(tsq(k)*qsq(k)) ), cov(k) )
1480 rprog_t(k,i,j,i_cov) = ( cov(k) * rho(k) - prog(k,i,j,i_cov) * dens(k,i,j) ) / dt
1484 rprog_t(k,i,j,i_cov) = - atmos_phy_bl_mynn_dump_coef * prog(k,i,j,i_cov) * dens(k,i,j) / dt
1489 nu(ks-1,i,j) = 0.0_rp
1490 kh(ks-1,i,j) = 0.0_rp
1495 prod(k,i,j) = 0.0_rp
1496 diss(k,i,j) = 0.0_rp
1501 dudz2(k,i,j) = undef
1504 cldfrac(k,i,j) = undef
1513 call file_history_query(hist_ri, do_put)
1514 if ( do_put )
call file_history_put(hist_ri, ri(:,:,:))
1515 call file_history_query(hist_pr, do_put)
1516 if ( do_put )
call file_history_put(hist_pr, pr(:,:,:))
1517 call file_history_query(hist_tke_pr, do_put)
1518 if ( do_put )
call file_history_put(hist_tke_pr, prod(:,:,:))
1519 call file_history_query(hist_tke_di, do_put)
1520 if ( do_put )
call file_history_put(hist_tke_di, diss(:,:,:))
1521 call file_history_query(hist_dudz2, do_put)
1522 if ( do_put )
call file_history_put(hist_dudz2, dudz2(:,:,:))
1523 call file_history_query(hist_lmix, do_put)
1524 if ( do_put )
call file_history_put(hist_lmix, l(:,:,:))
1526 call file_history_query(hist_flxu, do_put)
1527 if ( do_put )
call file_history_put(hist_flxu, flxu(1:,:,:))
1528 call file_history_query(hist_flxv, do_put)
1529 if ( do_put )
call file_history_put(hist_flxv, flxv(1:,:,:))
1530 call file_history_query(hist_flxt, do_put)
1531 if ( do_put )
call file_history_put(hist_flxt, flxt(1:,:,:))
1532 call file_history_query(hist_flxq, do_put)
1533 if ( do_put )
call file_history_put(hist_flxq, flxq(1:,:,:))
1535 call file_history_query(hist_flxu2, do_put)
1536 if ( do_put )
call file_history_put(hist_flxu2, flxu2(1:,:,:))
1537 call file_history_query(hist_flxv2, do_put)
1538 if ( do_put )
call file_history_put(hist_flxv2, flxv2(1:,:,:))
1539 call file_history_query(hist_flxt2, do_put)
1540 if ( do_put )
call file_history_put(hist_flxt2, flxt2(1:,:,:))
1541 call file_history_query(hist_flxq2, do_put)
1542 if ( do_put )
call file_history_put(hist_flxq2, flxq2(1:,:,:))
1553 subroutine mynn_main( &
1556 tke, tsq, qsq, cov, &
1558 Nu, RHONu, Kh, RHOKh, Pr, &
1559 prod, prod_t, prod_q, prod_c, &
1561 sm25, smp, sh25, shpgh, &
1563 dudz2, n2_new, Ri, &
1566 uflux, vflux, tflux, qflux, eflux, &
1567 Qlp, cldfrac, Zi, SFLX_BUOY, &
1573 SFLX_PT, SFLX_SH, SFLX_QV, SFC_DENS, &
1579 mynn_level3, initialize, &
1583 PTLV, Nu_f, Kh_f, q2_2, ac, betat, betaq, &
1585 f_smp, f_shpgh, f_gamma, tltv25, qwtv25, tvsq25, &
1586 tvsq_up, tvsq_lo, dtsq, dqsq, dcov, &
1588 Uh, Vh, Wh, qw2, qh, tlh, tvh, eh, dh, &
1597 matrix_solver_tridiagonal
1598 integer,
intent(in) :: KA, KS, KE_PBL
1599 integer,
intent(in) :: i, j
1600 real(RP),
intent(inout) :: tke(KA)
1601 real(RP),
intent(inout) :: tsq(KA)
1602 real(RP),
intent(inout) :: qsq(KA)
1603 real(RP),
intent(inout) :: cov(KA)
1604 real(RP),
intent(out) :: q(KA)
1605 real(RP),
intent(out) :: l(KA)
1606 real(RP),
intent(out) :: lq(KA)
1607 real(RP),
intent(out) :: Nu(KA)
1608 real(RP),
intent(out) :: RHONu(KA)
1609 real(RP),
intent(out) :: Kh(KA)
1610 real(RP),
intent(out) :: RHOKh(KA)
1611 real(RP),
intent(out) :: Pr(KA)
1612 real(RP),
intent(out) :: prod(KA)
1613 real(RP),
intent(out) :: prod_t(KA)
1614 real(RP),
intent(out) :: prod_q(KA)
1615 real(RP),
intent(out) :: prod_c(KA)
1616 real(RP),
intent(out) :: diss(KA)
1617 real(RP),
intent(out) :: diss_p(KA)
1618 real(RP),
intent(out) :: sm25(KA)
1619 real(RP),
intent(out) :: smp(KA)
1620 real(RP),
intent(out) :: sh25(KA)
1621 real(RP),
intent(out) :: shpgh(KA)
1622 real(RP),
intent(out) :: gammat(KA)
1623 real(RP),
intent(out) :: gammaq(KA)
1624 real(RP),
intent(out) :: dudz2(KA)
1625 real(RP),
intent(out) :: n2_new(KA)
1626 real(RP),
intent(out) :: Ri(KA)
1627 real(RP),
intent(out) :: dtldz(KA)
1628 real(RP),
intent(out) :: dqwdz(KA)
1629 real(RP),
intent(out) :: RHO(KA)
1630 real(RP),
intent(out) :: RHO_h(KA)
1631 real(RP),
intent(out) :: uflux(KA)
1632 real(RP),
intent(out) :: vflux(KA)
1633 real(RP),
intent(out) :: tflux(KA)
1634 real(RP),
intent(out) :: qflux(KA)
1635 real(RP),
intent(out) :: eflux(KA)
1636 real(RP),
intent(out) :: Qlp(KA)
1637 real(RP),
intent(out) :: cldfrac(KA)
1638 real(RP),
intent(out) :: Zi
1639 real(RP),
intent(out) :: SFLX_BUOY
1640 real(RP),
intent(in) :: U(KA)
1641 real(RP),
intent(in) :: V(KA)
1642 real(RP),
intent(in) :: W(KA)
1643 real(RP),
intent(in) :: DENS(KA)
1644 real(RP),
intent(in) :: PRES(KA)
1645 real(RP),
intent(in) :: POTT(KA)
1646 real(RP),
intent(in) :: POTL(KA)
1647 real(RP),
intent(in) :: POTV(KA)
1648 real(RP),
intent(in) :: Qw(KA)
1649 real(RP),
intent(in) :: N2(KA)
1650 real(RP),
intent(in) :: EXNER(KA)
1651 real(RP),
intent(in) :: QDRY(KA)
1652 real(RP),
intent(in) :: SFLX_PT
1653 real(RP),
intent(in) :: SFLX_SH
1654 real(RP),
intent(in) :: SFLX_QV
1655 real(RP),
intent(in) :: SFC_DENS
1656 real(RP),
intent(in) :: RLmo
1657 real(RP),
intent(in) :: us
1658 real(RP),
intent(in) :: ts
1659 real(RP),
intent(in) :: qs
1660 real(RP),
intent(in) :: z(KA)
1661 real(RP),
intent(in) :: CDZ(KA)
1662 real(RP),
intent(in) :: FDZ(KA)
1663 real(RP),
intent(in) :: F2H(KA,2)
1664 real(RP),
intent(in) :: frac_land
1665 real(RP),
intent(in) :: dt
1666 integer,
intent(in) :: I_B_TYPE
1667 logical,
intent(in) :: mynn_level3
1668 logical,
intent(in) :: initialize
1671 real(RP),
intent(out) :: PTLV (KA)
1672 real(RP),
intent(out) :: Nu_f (KA)
1673 real(RP),
intent(out) :: Kh_f (KA)
1674 real(RP),
intent(out) :: q2_2 (KA)
1675 real(RP),
intent(out) :: ac (KA)
1676 real(RP),
intent(out) :: betat(KA)
1677 real(RP),
intent(out) :: betaq(KA)
1679 real(RP),
intent(out) :: flx(0:KA)
1680 real(RP),
intent(out) :: a(KA)
1681 real(RP),
intent(out) :: b(KA)
1682 real(RP),
intent(out) :: c(KA)
1683 real(RP),
intent(out) :: d(KA)
1688 real(RP) :: phi_m, phi_h
1692 real(RP),
intent(out) :: f_smp (KA)
1693 real(RP),
intent(out) :: f_shpgh(KA)
1694 real(RP),
intent(out) :: f_gamma(KA)
1695 real(RP),
intent(out) :: tltv25 (KA)
1696 real(RP),
intent(out) :: qwtv25 (KA)
1697 real(RP),
intent(out) :: tvsq25 (KA)
1698 real(RP),
intent(out) :: tvsq_up(KA)
1699 real(RP),
intent(out) :: tvsq_lo(KA)
1700 real(RP),
intent(out) :: dtsq (KA)
1701 real(RP),
intent(out) :: dqsq (KA)
1702 real(RP),
intent(out) :: dcov (KA)
1710 real(RP),
intent(out) :: mflux(0:KA)
1713 real(RP),
intent(out) :: Uh(KA), Vh(KA), Wh(KA)
1714 real(RP),
intent(out) :: qw2(KA), qh(KA)
1715 real(RP),
intent(out) :: tlh(KA), tvh(KA)
1716 real(RP),
intent(out) :: eh(KA), dh(KA)
1717 real(RP),
intent(out) :: TEML(KA), LHVL(KA), psat(KA)
1720 real(RP),
intent(out) :: work(KA,4)
1725 integer :: k, it, nit
1728 if ( atmos_phy_bl_mynn_similarity .or. initialize )
then
1729 call get_phi( zeta, phi_m, phi_h, &
1730 z(ks), rlmo, i_b_type )
1733 call calc_vertical_differece( ka, ks, ke_pbl, &
1735 dudz2(:), dtldz(:), dqwdz(:), &
1736 u(:), v(:), potl(:), &
1738 cdz(:), fdz(:), f2h(:,:), &
1740 phi_m, phi_h, z(ks), &
1741 uh(:), vh(:), qw2(:), qh(:) )
1746 q(k) = sqrt( max( tke(k), atmos_phy_bl_mynn_tke_min ) * 2.0_rp )
1749 if ( atmos_phy_bl_mynn_use_zi )
then
1751 ptlv(k) = potl(k) * ( 1.0_rp + epstvap * qw(k) )
1753 call calc_zi_o2019( ka, ks, ke_pbl, &
1760 if ( initialize .or. (.not. mynn_level3) )
then
1764 n2_new(k) = min( max( n2(k), - atmos_phy_bl_mynn_n2_max ), atmos_phy_bl_mynn_n2_max )
1765 ri(k) = n2_new(k) / dudz2(k)
1767 if ( initialize )
then
1768 dudz2(ks) = max( dudz2(ks), 1e-4_rp )
1769 n2_new(ks) = min( n2_new(ks), 0.0_rp )
1770 ri(ks) = n2_new(ks) / dudz2(ks)
1773 sflx_buoy = - us3 * rlmo / karman
1775 if ( atmos_phy_bl_mynn_mf )
then
1781 dens(:), potv(:), potl(:), qw(:), &
1782 u(:), v(:), w(:), tke(:), &
1785 sflx_sh, sflx_buoy, &
1788 zi, z(:), cdz(:), f2h(:,:), &
1790 tflux(:), qflux(:), &
1791 uflux(:), vflux(:), eflux(:), &
1792 tlh(:), tvh(:), qh(:), &
1793 uh(:), vh(:), wh(:), eh(:), dh(:) )
1804 z(:), cdz(:), fdz(:), &
1808 call get_q2_level2( &
1810 dudz2(:), ri(:), l(:), &
1813 if ( initialize )
then
1815 q(k) = sqrt( q2_2(k) )
1820 ac(k) = min( q(k) / sqrt( q2_2(k) + 1e-20_rp ), 1.0_rp )
1829 potv(:), dudz2(:), &
1830 dtldz(:), dqwdz(:), &
1831 betat(:), betaq(:), &
1833 tsq(:), qsq(:), cov(:), &
1834 sm25(:), f_smp(:), &
1835 sh25(:), f_shpgh(:), &
1837 tltv25(:), qwtv25(:), &
1839 tvsq_up(:), tvsq_lo(:) )
1844 flx(ke_pbl) = 0.0_rp
1846 if ( initialize )
then
1855 call partial_condensation( ka, ks, ke_pbl, &
1856 betat(:), betaq(:), &
1857 qlp(:), cldfrac(:), &
1861 tsq(:), qsq(:), cov(:), &
1862 teml(:), lhvl(:), psat(:) )
1866 n2_new(k) = min(atmos_phy_bl_mynn_n2_max, &
1867 grav * ( dtldz(k) * betat(k) + dqwdz(k) * betaq(k) ) / potv(k) )
1868 ri(k) = n2_new(k) / dudz2(k)
1871 sflx_buoy = grav / potv(ks) * ( betat(ks) * sflx_pt + betaq(ks) * sflx_qv ) / sfc_dens
1874 if ( atmos_phy_bl_mynn_use_zi )
then
1876 ptlv(k) = potl(k) * betat(k) + qw(k) * betaq(k)
1878 call calc_zi_o2019( ka, ks, ke_pbl, &
1885 if ( atmos_phy_bl_mynn_mf )
then
1888 dens(:), potv(:), potl(:), qw(:), &
1889 u(:), v(:), w(:), tke(:), &
1892 sflx_sh, sflx_buoy, &
1895 zi, z(:), cdz(:), f2h(:,:), &
1897 tflux(:), qflux(:), &
1898 uflux(:), vflux(:), eflux(:), &
1899 tlh(:), tvh(:), qh(:), &
1900 uh(:), vh(:), wh(:), eh(:), dh(:) )
1912 z(:), cdz(:), fdz(:), &
1916 call get_q2_level2( &
1918 dudz2(:), ri(:), l(:), &
1922 ac(k) = min( q(k) / sqrt( q2_2(k) + 1e-20_rp ), 1.0_rp )
1931 potv(:), dudz2(:), &
1932 dtldz(:), dqwdz(:), &
1933 betat(:), betaq(:), &
1935 initialize .and. it==1, &
1936 tsq(:), qsq(:), cov(:), &
1937 sm25(:), f_smp(:), &
1938 sh25(:), f_shpgh(:), &
1940 tltv25(:), qwtv25(:), &
1942 tvsq_up(:), tvsq_lo(:) )
1946 nu_f(k) = lq(k) * sm25(k)
1947 kh_f(k) = lq(k) * sh25(k)
1949 if ( atmos_phy_bl_mynn_similarity )
then
1950 nu_f(ks) = karman * z(ks) * us / phi_m
1951 kh_f(ks) = karman * z(ks) * us / phi_h
1955 nu(k) = min( f2h(k,1) * nu_f(k+1) + f2h(k,2) * nu_f(k), &
1956 atmos_phy_bl_mynn_nu_max )
1957 kh(k) = min( f2h(k,1) * kh_f(k+1) + f2h(k,2) * kh_f(k), &
1958 atmos_phy_bl_mynn_kh_max )
1962 sw = 0.5_rp - sign(0.5_rp, abs(kh(k)) - eps)
1963 pr(k) = nu(k) / ( kh(k) + sw ) * ( 1.0_rp - sw ) &
1967 rho(ks) = dens(ks) + dt * sflx_qv / cdz(ks)
1973 rho_h(k) = f2h(k,1) * rho(k+1) + f2h(k,2) * rho(k)
1974 rhonu(k) = max( nu(k) * rho_h(k), 1e-20_rp )
1975 rhokh(k) = max( kh(k) * rho_h(k), 1e-20_rp )
1981 if ( mynn_level3 )
then
1985 tltv = betat(k) * tsq(k) + betaq(k) * cov(k)
1986 qwtv = betat(k) * cov(k) + betaq(k) * qsq(k)
1987 gammat(k) = f_gamma(k) * ( tltv - tltv25(k) )
1988 gammaq(k) = f_gamma(k) * ( qwtv - qwtv25(k) )
1990 wtl = - lq(k) * ( sh25(k) * dtldz(k) + gammat(k) )
1991 wqw = - lq(k) * ( sh25(k) * dqwdz(k) + gammaq(k) )
1992 prod_t(k) = - 2.0_rp * wtl * dtldz(k)
1993 prod_q(k) = - 2.0_rp * wqw * dqwdz(k)
1994 prod_c(k) = - wtl * dqwdz(k) - wqw * dtldz(k)
1997 if ( atmos_phy_bl_mynn_similarity .or. initialize )
then
1998 tmp = 2.0_rp * us * phi_h / ( karman * z(ks) )
1999 tmp = tmp * ( zeta / ( z(ks) * rlmo ) )**2
2001 prod_t(ks) = tmp * ts**2
2003 prod_q(ks) = tmp * ts * qs
2005 prod_c(ks) = tmp * qs**2
2010 flx(ke_pbl) = 0.0_rp
2012 flx(k) = rhonu(k) * ( tsq(k+1) - tsq(k) ) / fdz(k)
2015 prod_t(k) = prod_t(k) + ( flx(k) - flx(k-1) ) / cdz(k)
2018 flx(k) = rhonu(k) * ( qsq(k+1) - qsq(k) ) / fdz(k)
2021 prod_q(k) = prod_q(k) + ( flx(k) - flx(k-1) ) / cdz(k)
2024 flx(k) = rhonu(k) * ( cov(k+1) - cov(k) ) / fdz(k)
2027 prod_c(k) = prod_c(k) + ( flx(k) - flx(k-1) ) / cdz(k)
2030 call get_gamma_implicit( &
2033 tsq(:), qsq(:), cov(:), &
2034 dtldz(:), dqwdz(:), potv(:), &
2035 prod_t(:), prod_q(:), prod_c(:), &
2036 betat(:), betaq(:), &
2037 f_gamma(:), l(:), q(:), &
2039 dtsq(:), dqsq(:), dcov(:) )
2043 tltv = betat(k) * ( tsq(k) + dtsq(k) ) + betaq(k) * ( cov(k) + dcov(k) )
2044 qwtv = betat(k) * ( cov(k) + dcov(k) ) + betaq(k) * ( qsq(k) + dqsq(k) )
2045 gammat(k) = f_gamma(k) * ( tltv - tltv25(k) )
2046 gammaq(k) = f_gamma(k) * ( qwtv - qwtv25(k) )
2048 tvsq = max( betat(k) * tltv + betaq(k) * qwtv, 0.0_rp )
2049 tvsq = tvsq - tvsq25(k)
2050 tvsq = min( max( tvsq, tvsq_lo(k) ), tvsq_up(k) )
2051 smp(k) = f_smp(k) * tvsq
2052 shpgh(k) = f_shpgh(k) * tvsq
2054 wtl = - lq(k) * ( sh25(k) * dtldz(k) + gammat(k) )
2055 wqw = - lq(k) * ( sh25(k) * dqwdz(k) + gammaq(k) )
2056 prod_t(k) = - 2.0_rp * wtl * dtldz(k)
2057 prod_q(k) = - 2.0_rp * wqw * dqwdz(k)
2058 prod_c(k) = - wtl * dqwdz(k) - wqw * dtldz(k)
2061 if ( atmos_phy_bl_mynn_similarity .or. initialize )
then
2067 tmp = 2.0_rp * us * phi_h / ( karman * z(ks) )
2068 tmp = tmp * ( zeta / ( z(ks) * rlmo ) )**2
2070 prod_t(ks) = tmp * ts**2
2072 prod_q(ks) = tmp * ts * qs
2074 prod_c(ks) = tmp * qs**2
2091 prod(k) = lq(k) * ( ( sm25(k) + smp(k) ) * dudz2(k) &
2092 - ( sh25(k) * n2_new(k) - shpgh(k) ) )
2094 if ( atmos_phy_bl_mynn_similarity .or. initialize )
then
2095 prod(ks) = us3 / ( karman * z(ks) ) * ( phi_m - zeta )
2099 diss(k) = - 2.0_rp * q(k) / ( b1 * l(k) )
2101 diss_p(k) = dt * 2.0_rp * q(k) / ( b2 * l(k) )
2105 if ( .not. initialize )
exit
2120 d(k) = dt * ( - ( flx(k) - flx(k-1) ) / ( cdz(k) * rho(k) ) + prod(k) )
2126 ap = - dt * atmos_phy_bl_mynn_sq_fact * rhonu(k) / fdz(k)
2127 a(k) = ap / ( rho(k) * cdz(k) )
2130 b(k) = - a(k) - diss(k) * dt
2132 b(k) = - a(k) + dt * atmos_phy_bl_mynn_sq_fact * rhonu(k-1) / ( fdz(k-1) * rho(k) * cdz(k) ) - diss(k) * dt
2135 b(k) = - a(k) - c(k) - diss(k) * dt
2137 c(k+1) = ap / ( rho(k+1) * cdz(k+1) )
2140 b(ke_pbl) = - c(ke_pbl) - diss(ke_pbl) * dt
2142 call matrix_solver_tridiagonal( &
2147 a(:), b(:), c(:), d(:), &
2151 tke(k) = min( max( tke(k), atmos_phy_bl_mynn_tke_min ), 100.0_rp )
2153 tke(ke_pbl) = 0.0_rp
2157 q(k) = ( q(k) + sqrt( tke(k) * 2.0_rp ) ) * 0.5_rp
2161 if ( .not. mynn_level3 ) cycle
2166 d(k) = dt * prod_t(k)
2171 ap = - dt * rhonu(k) / fdz(k)
2172 a(k) = ap / ( rho(k) * cdz(k) )
2175 b(k) = - a(k) + diss_p(k)
2177 b(k) = - a(k) + dt * rhonu(k-1) / ( fdz(k-1) * rho(k) * cdz(k) ) + diss_p(k)
2180 b(k) = - a(k) - c(k) + 1.0_rp + diss_p(k)
2182 c(k+1) = ap / ( rho(k+1) * cdz(k+1) )
2185 b(ke_pbl) = - c(ke_pbl) + diss_p(ke_pbl)
2187 call matrix_solver_tridiagonal( &
2192 a(:), b(:), c(:), d(:), &
2196 tsq(k) = max( tsq(k), 0.0_rp )
2198 tsq(ke_pbl) = 0.0_rp
2204 d(k) = dt * prod_q(k)
2209 call matrix_solver_tridiagonal( &
2214 a(:), b(:), c(:), d(:), &
2218 qsq(k) = max( qsq(k), 0.0_rp )
2220 qsq(ke_pbl) = 0.0_rp
2226 d(k) = dt * prod_c(k)
2231 call matrix_solver_tridiagonal( &
2236 a(:), b(:), c(:), d(:), &
2240 cov(k) = sign( min( abs(cov(k)), sqrt(tsq(k)*qsq(k)) ), cov(k) )
2242 cov(ke_pbl) = 0.0_rp
2247 end subroutine mynn_main
2250 subroutine get_length( &
2266 integer,
intent(in),
value :: KA, KS, KE_PBL
2268 real(RP),
intent(in) :: q(KA)
2269 real(RP),
intent(in) :: n2(KA)
2270 real(RP),
intent(in) :: POTV(KA)
2271 real(RP),
intent(in) :: DENS(KA)
2272 real(RP),
intent(in) :: mflux_gl(0:KA)
2273 real(RP),
intent(in) :: Zi
2274 real(RP),
intent(in) :: SFLX_BUOY
2275 real(RP),
intent(in) :: RLmo
2276 real(RP),
intent(in) :: z(KA)
2277 real(RP),
intent(in) :: CDZ(KA)
2278 real(RP),
intent(in) :: FDZ(KA)
2279 logical,
intent(in) :: initialize
2281 real(RP),
intent(out) :: l(KA)
2283 real(RP),
parameter :: sqrt05 = sqrt( 0.5_rp )
2306 if ( atmos_phy_bl_mynn_use_zi )
then
2310 if ( z(k) > zi * 1.3_rp )
then
2315 kmax = max(kmax, ks)
2325 int_qz = int_qz + z(k) * qdz
2329 lt = min( max(0.23_rp * int_qz / (int_q + 1e-20_rp), &
2331 atmos_phy_bl_mynn_lt_max )
2334 qc = ( lt * max(sflx_buoy,0.0_rp) )**oneoverthree
2335 if ( atmos_phy_bl_mynn_o2019 )
then
2336 tau = 0.5_rp * zi / max(sflx_buoy,1.0e-10_rp)**oneoverthree
2344 sw = sign(0.5_rp, zeta) + 0.5_rp
2345 ls = karman * z(k) &
2346 * ( sw / (1.0_rp + atmos_phy_bl_mynn_cns*zeta*sw ) &
2347 + ( (1.0_rp - atmos_phy_bl_mynn_alpha4*zeta)*(1.0_rp-sw) )**0.2_rp )
2350 sw = sign(0.5_rp, n2(k)-eps) + 0.5_rp
2351 rn2sr = 1.0_rp / ( sqrt(n2(k)*sw) + 1.0_rp-sw)
2352 if ( atmos_phy_bl_mynn_o2019 )
then
2353 if ( z(k) > zi ) tau = 50.0_rp
2375 lb = atmos_phy_bl_mynn_alpha2 * max( q(k), ( mflux_gl(k)+mflux_gl(k-1))/dens(k) ) * rn2sr * sw &
2376 + tau * q(k) * sqrt05 * (1.0_rp-sw)
2378 lb = atmos_phy_bl_mynn_alpha2 * (1.0_rp + 5.0_rp * sqrt(qc*rn2sr*rlt)) * q(k) * rn2sr * sw &
2379 + 1.e10_rp * (1.0_rp-sw)
2383 if ( initialize )
then
2384 l(k) = min( ls, lb )
2385 else if ( atmos_phy_bl_mynn_o2019 )
then
2386 l(k) = min( 1.0_rp / ( 1.0_rp/ls + rlt ), lb )
2388 l(k) = 1.0_rp / ( 1.0_rp/ls + rlt + 1.0_rp/(lb+1e-20_rp) )
2393 end subroutine get_length
2396 subroutine get_q2_level2( &
2402 integer,
intent(in),
value :: KA, KS, KE_PBL
2404 real(RP),
intent(in) :: dudz2(KA)
2405 real(RP),
intent(in) :: Ri(KA)
2406 real(RP),
intent(in) :: l(KA)
2408 real(RP),
intent(out) :: q2_2(KA)
2415 real(RP) :: A2k, F1, Rf1, AF12
2420 if ( atmos_phy_bl_mynn_k2010 .and. ri(k) > 0.0_rp )
then
2421 a2k = a2 / ( 1.0_rp + ri(k) )
2425 f1 = b1 * ( g1 - c1 ) + 2.0 * a1 * ( 3.0_rp - 2.0_rp * c2 ) + 3.0 * a2k * ( 1.0_rp - c2 ) * ( 1.0_rp - c5 )
2426 rf1 = b1 * ( g1 - c1 ) / f1
2427 af12 = a1 * f1 / ( a2k * f2 )
2428 rf = min(0.5_rp / af12 * ( ri(k) &
2430 - sqrt(ri(k)**2 + 2.0_rp*af12*(rf1-2.0_rp*rf2)*ri(k) + (af12*rf1)**2) ), &
2432 sh_2 = 3.0_rp * a2k * (g1+g2) * (rfc-rf) / (1.0_rp-rf)
2433 sm_2 = sh_2 * af12 * (rf1-rf) / (rf2-rf)
2434 q2_2(k) = b1 * l(k)**2 * sm_2 * (1.0_rp-rf) * dudz2(k)
2438 end subroutine get_q2_level2
2441 subroutine get_smsh( &
2464 integer,
intent(in),
value :: KA, KS, KE_PBL
2465 integer,
intent(in),
value :: i, j
2467 real(RP),
intent(in) :: q(KA)
2468 real(RP),
intent(in) :: ac(KA)
2469 real(RP),
intent(in) :: l(KA)
2470 real(RP),
intent(in) :: n2(KA)
2471 real(RP),
intent(in) :: Ri(KA)
2472 real(RP),
intent(in) :: potv(KA)
2473 real(RP),
intent(in) :: dudz2(KA)
2474 real(RP),
intent(in) :: dtldz(KA)
2475 real(RP),
intent(in) :: dqwdz(KA)
2476 real(RP),
intent(in) :: betat(KA)
2477 real(RP),
intent(in) :: betaq(KA)
2478 logical,
intent(in),
value :: mynn_level3
2479 logical,
intent(in),
value :: initialize
2481 real(RP),
intent(inout) :: tsq(KA)
2482 real(RP),
intent(inout) :: qsq(KA)
2483 real(RP),
intent(inout) :: cov(KA)
2485 real(RP),
intent(out) :: sm25 (KA)
2486 real(RP),
intent(out) :: f_smp (KA)
2487 real(RP),
intent(out) :: sh25 (KA)
2488 real(RP),
intent(out) :: f_shpgh(KA)
2489 real(RP),
intent(out) :: f_gamma(KA)
2490 real(RP),
intent(out) :: tltv25 (KA)
2491 real(RP),
intent(out) :: qwtv25 (KA)
2492 real(RP),
intent(out) :: tvsq25 (KA)
2493 real(RP),
intent(out) :: tvsq_up(KA)
2494 real(RP),
intent(out) :: tvsq_lo(KA)
2508 real(RP) :: f1, f2, f3, f4
2523 real(RP) :: tmp1, tmp2
2533 if ( atmos_phy_bl_mynn_k2010 .and. ri(k) > 0.0_rp )
then
2534 a2k = a2 / ( 1.0_rp + ri(k) )
2543 f1 = - 3.0_rp * ac2 * a2k * b2 * ( 1.0_rp - c3 )
2544 f2 = - 9.0_rp * ac2 * a1 * a2k * ( 1.0_rp - c2 )
2545 f3 = 9.0_rp * ac2 * a2k**2 * ( 1.0_rp - c2 ) * ( 1.0_rp - c5 )
2546 f4 = - 12.0_rp * ac2 * a1 * a2k * ( 1.0_rp - c2 )
2548 if ( mynn_level3 )
then
2549 ghq2 = max( - n2(k) * l2, -q2 )
2553 gmq2 = dudz2(k) * l2
2555 p1q2 = q2 + f1 * ghq2
2556 p2q2 = q2 + f2 * ghq2
2557 p3q2 = q2 + ( f1 + f3 ) * ghq2
2558 p4q2 = q2 + ( f1 + f4 ) * ghq2
2559 p5q2 = 6.0_rp * ac2 * a1**2 * gmq2
2560 rd25q2 = q2 / max( p2q2 * p4q2 + p5q2 * p3q2, 1e-20_rp )
2562 sm25(k) = ac(k) * a1 * ( p3q2 - 3.0_rp * c1 * p4q2 ) * rd25q2
2563 sh25(k) = ac(k) * a2k * ( p2q2 + 3.0_rp * c1 * p5q2 ) * rd25q2
2565 fact = ac(k) * b2 * l2 * sh25(k)
2566 tsq25 = fact * dtldz(k)**2
2567 qsq25 = fact * dqwdz(k)**2
2568 cov25 = fact * dtldz(k) * dqwdz(k)
2570 if ( initialize .or. (.not. mynn_level3 ) )
then
2576 if ( mynn_level3 )
then
2578 if ( q2 <= 1e-10_rp )
then
2588 tltv25(k) = betat(k) * tsq25 + betaq(k) * cov25
2589 qwtv25(k) = betat(k) * cov25 + betaq(k) * qsq25
2590 tvsq25(k) = max(betat(k) * tltv25(k) + betaq(k) * qwtv25(k), 0.0_rp)
2591 cw25 = p1q2 * ( p2q2 + 3.0_rp * c1 * p5q2 ) * rd25q2 / ( 3.0_rp * q2 )
2594 l2q2 = min( l2q2, 1.0_rp / max(n2(k), eps) )
2597 rdpq2 = q2 / max( p2q2 * ( f4 * ghq2 + q2 ) + p5q2 * ( f3 * ghq2 + q2 ), 1e-20_rp )
2599 tltv = betat(k) * tsq(k) + betaq(k) * cov(k)
2600 qwtv = betat(k) * cov(k) + betaq(k) * qsq(k)
2601 tvsq = max(betat(k) * tltv + betaq(k) * qwtv, 0.0_rp)
2603 ew = ( 1.0_rp - c3 ) * ( - p2q2 * f4 - p5q2 * f3 ) * rdpq2
2604 if ( abs(ew) > eps )
then
2605 fact = q2 * potv(k)**2 / ( ew * l2q2 * grav**2 )
2606 if ( fact > 0.0_rp )
then
2613 tvsq_up(k) = ( tmp1 - cw25 ) * fact
2614 tvsq_lo(k) = ( tmp2 - cw25 ) * fact
2620 emq2 = 3.0_rp * ac(k) * a1 * ( 1.0_rp - c3 ) * ( f3 - f4 ) * rdpq2
2621 eh = 3.0_rp * ac(k) * a2k * ( 1.0_rp - c3 ) * ( p2q2 + p5q2 ) * rdpq2
2624 fact = grav / potv(k)
2625 f_smp(k) = emq2 * fact**2 * l2q2
2626 f_shpgh(k) = eh * fact**2 / q2
2627 f_gamma(k) = - eh * fact / q2
2640 end subroutine get_smsh
2643 subroutine get_gamma_implicit( &
2647 dtldz, dqwdz, POTV, &
2648 prod_t, prod_q, prod_c, &
2657 integer,
intent(in),
value :: KA, KS, KE
2658 integer,
intent(in),
value :: i, j
2660 real(RP),
intent(in) :: tsq (KA)
2661 real(RP),
intent(in) :: qsq (KA)
2662 real(RP),
intent(in) :: cov (KA)
2663 real(RP),
intent(in) :: dtldz (KA)
2664 real(RP),
intent(in) :: dqwdz (KA)
2665 real(RP),
intent(in) :: POTV (KA)
2666 real(RP),
intent(in) :: prod_t (KA)
2667 real(RP),
intent(in) :: prod_q (KA)
2668 real(RP),
intent(in) :: prod_c (KA)
2669 real(RP),
intent(in) :: betat (KA)
2670 real(RP),
intent(in) :: betaq (KA)
2671 real(RP),
intent(in) :: f_gamma(KA)
2672 real(RP),
intent(in) :: l (KA)
2673 real(RP),
intent(in) :: q (KA)
2674 real(RP),
intent(in),
value :: dt
2676 real(RP),
intent(out) :: dtsq(KA)
2677 real(RP),
intent(out) :: dqsq(KA)
2678 real(RP),
intent(out) :: dcov(KA)
2680 real(RP) :: a11, a12, a21, a22, a23, a32, a33
2681 real(RP) :: v1, v2, v3
2683 real(RP) :: a13, a31, det
2692 f1 = f_gamma(k) * l(k) * q(k)
2693 f2 = - 2.0_rp * q(k) / ( b2 * l(k) )
2695 v1 = prod_t(k) + f2 * tsq(k)
2696 v2 = prod_c(k) + f2 * cov(k)
2697 v3 = prod_q(k) + f2 * qsq(k)
2699 a11 = - f1 * dtldz(k) * betat(k) * 2.0_rp + 1.0_rp / dt - f2
2700 a12 = - f1 * dtldz(k) * betaq(k) * 2.0_rp
2701 a21 = - f1 * dqwdz(k) * betat(k)
2702 a22 = - f1 * ( dtldz(k) * betat(k) + dqwdz(k) * betaq(k) ) + 1.0_rp / dt - f2
2703 a23 = - f1 * dtldz(k) * betaq(k)
2704 a32 = - f1 * dqwdz(k) * betat(k) * 2.0_rp
2705 a33 = - f1 * dqwdz(k) * betaq(k) * 2.0_rp + 1.0_rp / dt - f2
2707 if ( q(k) < 1e-20_rp .or. abs(f_gamma(k)) * dt >= q(k) )
then
2711 dcov(k) = ( v2 - f1 * v1 - f2 * v3 ) / ( a22 - f1 * a12 - f2 * a32 )
2712 dtsq(k) = ( v1 - a12 * dcov(k) ) / a11
2713 dqsq(k) = ( v3 - a32 * dcov(k) ) / a33
2717 f1 = - l(k) * grav / potv(k) * f_gamma(k) / q(k)
2718 f2 = f1 * betat(k)**2
2722 f2 = f1 * betat(k) * betaq(k) * 2.0_rp
2726 f2 = f1 * betaq(k)**2
2731 det = a11 * ( a22 * a33 - a23 * a32 ) + a12 * ( a23 * a31 - a21 * a33 ) + a13 * ( a21 * a32 - a22 * a31 )
2732 dtsq(k) = ( v1 * ( a22 * a33 - a23 * a32 ) + a12 * ( a23 * v3 - v2 * a33 ) + a13 * ( v2 * a32 - a22 * v3 ) ) / det
2733 dcov(k) = ( a11 * ( v2 * a33 - a23 * v3 ) + v1 * ( a23 * a31 - a21 * a33 ) + a13 * ( a21 * v3 - v2 * a31 ) ) / det
2734 dqsq(k) = ( a11 * ( a22 * v3 - v2 * a32 ) + a12 * ( v2 * a31 - a21 * v3 ) + v3 * ( a21 * a32 - a22 * a31 ) ) / det