--- trunk/phylmd/coef_diff_turb.f 2018/02/05 10:39:38 254 +++ trunk/phylmd/Interface_surf/coef_diff_turb.f 2018/09/06 13:19:51 302 @@ -4,23 +4,23 @@ contains - subroutine coef_diff_turb(dtime, nsrf, ni, paprs, pplay, u, v, q, t, ts, & - cdragm, zgeop, coefm, coefh, q2) + subroutine coef_diff_turb(nsrf, ni, paprs, pplay, u, v, q, t, ts, cdragm, & + zgeop, coefm, coefh, q2) ! Computes coefficients for turbulent diffusion in the atmosphere. + use nr_util, only: assert + USE clesphys, ONLY: ok_kzmin use coefkz_m, only: coefkz use coefkzmin_m, only: coefkzmin use coefkz2_m, only: coefkz2 USE conf_phys_m, ONLY: iflag_pbl USE dimphy, ONLY: klev - use nr_util, only: assert USE suphec_m, ONLY: rd, rg, rkappa use ustarhb_m, only: ustarhb use yamada4_m, only: yamada4 - REAL, INTENT(IN):: dtime ! interval du temps (secondes) INTEGER, INTENT(IN):: nsrf INTEGER, INTENT(IN):: ni(:) ! (knon) REAL, INTENT(IN):: paprs(:, :) ! (knon, klev + 1) @@ -39,7 +39,6 @@ REAL coefm0(size(ni), 2:klev), coefh0(size(ni), 2:klev) ! (knon, 2:klev) REAL zlay(size(ni), klev), teta(size(ni), klev) ! (knon, klev) real zlev(size(ni), klev + 1) ! (knon, klev + 1) - REAL ustar(size(ni)) ! (knon) integer k !------------------------------------------------------------------------- @@ -48,22 +47,6 @@ size(v, 1), size(q, 1), size(t, 1), size(ts), size(cdragm), & size(zgeop, 1), size(coefm, 1), size(coefh, 1), size(q2, 1)], & "coef_diff_turb knon") - - CALL coefkz(nsrf, paprs, pplay, ts, u, v, t, q, zgeop, coefm, coefh) - - IF (iflag_pbl == 1) THEN - CALL coefkz2(nsrf, paprs, pplay, t, coefm0, coefh0) - coefm = max(coefm, coefm0) - coefh = max(coefh, coefh0) - END IF - - IF (ok_kzmin) THEN - ! Calcul d'une diffusion minimale pour les conditions tres stables - CALL coefkzmin(paprs, pplay, u, v, t, q, cdragm, coefh0) - coefm0 = coefh0 - coefm = max(coefm, coefm0) - coefh = max(coefh, coefh0) - END IF IF (iflag_pbl >= 6) THEN ! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et @@ -76,20 +59,30 @@ / paprs(:, k) * (pplay(:, k-1) - pplay(:, k)) / rg END DO - DO k = 1, klev - teta(:, k) = t(:, k) * (paprs(:, 1) / pplay(:, k))**rkappa & - * (1. + 0.61 * q(:, k)) - END DO + forall (k = 1:klev) teta(:, k) = t(:, k) & + * (paprs(:, 1) / pplay(:, k))**rkappa * (1. + 0.61 * q(:, k)) zlev(:, 1) = 0. + forall (k = 2:klev) zlev(:, k) = 0.5 * (zlay(:, k) + zlay(:, k-1)) zlev(:, klev + 1) = 2. * zlay(:, klev) - zlay(:, klev - 1) - DO k = 2, klev - zlev(:, k) = 0.5 * (zlay(:, k) + zlay(:, k-1)) - END DO - - ustar = ustarhb(u(:, 1), v(:, 1), cdragm) - CALL yamada4(dtime, zlev, zlay, u, v, teta, q2, coefm, coefh, ustar) + CALL yamada4(zlev, zlay, u, v, teta, q2, coefm, coefh, & + ustarhb(u(:, 1), v(:, 1), cdragm)) + else + CALL coefkz(nsrf, paprs, pplay, ts, u, v, t, q, zgeop, coefm, coefh) + + IF (iflag_pbl == 1) THEN + CALL coefkz2(nsrf, paprs, pplay, t, coefm0, coefh0) + coefm = max(coefm, coefm0) + coefh = max(coefh, coefh0) + END IF + + IF (ok_kzmin) THEN + ! Calcul d'une diffusion minimale pour les conditions tres stables + CALL coefkzmin(paprs, pplay, u, v, t, q, cdragm, coefh0) + coefm = max(coefm, coefh0) + coefh = max(coefh, coefh0) + END IF END IF end subroutine coef_diff_turb