/[lmdze]/trunk/phylmd/Interface_surf/coef_diff_turb.f
ViewVC logotype

Diff of /trunk/phylmd/Interface_surf/coef_diff_turb.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

trunk/Sources/phylmd/coef_diff_turb.f revision 250 by guez, Fri Jan 5 18:18:53 2018 UTC trunk/phylmd/Interface_surf/coef_diff_turb.f revision 288 by guez, Tue Jul 24 16:27:12 2018 UTC
# Line 4  module coef_diff_turb_m Line 4  module coef_diff_turb_m
4    
5  contains  contains
6    
7    subroutine coef_diff_turb(dtime, nsrf, ni, ypaprs, ypplay, yu, yv, yq, yt, &    subroutine coef_diff_turb(dtime, nsrf, ni, paprs, pplay, u, v, q, t, ts, &
8         yts, ycdragm, zgeop, ycoefm, ycoefh, yq2)         cdragm, zgeop, coefm, coefh, q2)
9    
10        ! Computes coefficients for turbulent diffusion in the atmosphere.
11    
12        use nr_util, only: assert
13    
14      USE clesphys, ONLY: ok_kzmin      USE clesphys, ONLY: ok_kzmin
15      use coefkz_m, only: coefkz      use coefkz_m, only: coefkz
16      use coefkzmin_m, only: coefkzmin      use coefkzmin_m, only: coefkzmin
17      use coefkz2_m, only: coefkz2      use coefkz2_m, only: coefkz2
18      USE conf_phys_m, ONLY: iflag_pbl      USE conf_phys_m, ONLY: iflag_pbl
19      USE dimphy, ONLY: klev, klon      USE dimphy, ONLY: klev
20      USE suphec_m, ONLY: rd, rg, rkappa      USE suphec_m, ONLY: rd, rg, rkappa
21      use ustarhb_m, only: ustarhb      use ustarhb_m, only: ustarhb
22      use yamada4_m, only: yamada4      use yamada4_m, only: yamada4
# Line 20  contains Line 24  contains
24      REAL, INTENT(IN):: dtime ! interval du temps (secondes)      REAL, INTENT(IN):: dtime ! interval du temps (secondes)
25      INTEGER, INTENT(IN):: nsrf      INTEGER, INTENT(IN):: nsrf
26      INTEGER, INTENT(IN):: ni(:) ! (knon)      INTEGER, INTENT(IN):: ni(:) ! (knon)
27      REAL, INTENT(IN):: ypaprs(:, :) ! (klon, klev + 1)      REAL, INTENT(IN):: paprs(:, :) ! (knon, klev + 1)
28      REAL, INTENT(IN):: ypplay(:, :) ! (klon, klev)      REAL, INTENT(IN):: pplay(:, :) ! (knon, klev)
29      REAL, INTENT(IN), dimension(:, :):: yu, yv, yq, yt ! (klon, klev)      REAL, INTENT(IN), dimension(:, :):: u, v, q, t ! (knon, klev)
30      REAL, INTENT(IN):: yts(:), ycdragm(:) ! (knon)      REAL, INTENT(IN):: ts(:), cdragm(:) ! (knon)
31      REAL, INTENT(IN):: zgeop(:, :) ! (knon, klev)      REAL, INTENT(IN):: zgeop(:, :) ! (knon, klev)
32      REAL, intent(out):: ycoefm(:, 2:) ! (knon, 2:klev) coefficient, vitesse      REAL, intent(out):: coefm(:, 2:) ! (knon, 2:klev) coefficient, vitesse
33    
34      real, intent(out):: ycoefh(:, 2:) ! (knon, 2:klev)      real, intent(out):: coefh(:, 2:) ! (knon, 2:klev)
35      ! coefficient, chaleur et humidité      ! coefficient, chaleur et humidité
36    
37      real, intent(inout):: yq2(:, :) ! (klon, klev + 1)      real, intent(inout):: q2(:, :) ! (knon, klev + 1)
38    
39      ! Local:      ! Local:
40      REAL ycoefm0(klon, 2:klev), ycoefh0(klon, 2:klev)      REAL coefm0(size(ni), 2:klev), coefh0(size(ni), 2:klev) ! (knon, 2:klev)
41      REAL yzlay(klon, klev), zlev(klon, klev + 1), yteta(klon, klev)      REAL zlay(size(ni), klev), teta(size(ni), klev) ! (knon, klev)
42      REAL ustar(klon)      real zlev(size(ni), klev + 1) ! (knon, klev + 1)
43      integer k, knon      integer k
44    
45      !-------------------------------------------------------------------------      !-------------------------------------------------------------------------
46    
47      knon = size(ni)      call assert(size(ni) == [size(paprs, 1), size(pplay, 1), size(u, 1), &
48      CALL coefkz(nsrf, ypaprs(:knon, :), ypplay(:knon, :), yts(:knon), &           size(v, 1), size(q, 1), size(t, 1), size(ts), size(cdragm), &
49           yu(:knon, :), yv(:knon, :), yt(:knon, :), yq(:knon, :), zgeop, &           size(zgeop, 1), size(coefm, 1), size(coefh, 1), size(q2, 1)], &
50           ycoefm, ycoefh)           "coef_diff_turb knon")
51        
52        CALL coefkz(nsrf, paprs, pplay, ts, u, v, t, q, zgeop, coefm, coefh)
53    
54      IF (iflag_pbl == 1) THEN      IF (iflag_pbl == 1) THEN
55         CALL coefkz2(nsrf, knon, ypaprs, ypplay, yt, ycoefm0(:knon, :), &         CALL coefkz2(nsrf, paprs, pplay, t, coefm0, coefh0)
56              ycoefh0(:knon, :))         coefm = max(coefm, coefm0)
57         ycoefm = max(ycoefm, ycoefm0(:knon, :))         coefh = max(coefh, coefh0)
        ycoefh = max(ycoefh, ycoefh0(:knon, :))  
58      END IF      END IF
59    
60      IF (ok_kzmin) THEN      IF (ok_kzmin) THEN
61         ! Calcul d'une diffusion minimale pour les conditions tres stables         ! Calcul d'une diffusion minimale pour les conditions tres stables
62         CALL coefkzmin(knon, ypaprs, ypplay, yu, yv, yt, yq, ycdragm(:knon), &         CALL coefkzmin(paprs, pplay, u, v, t, q, cdragm, coefh0)
63              ycoefh0(:knon, :))         coefm = max(coefm, coefh0)
64         ycoefm0(:knon, :) = ycoefh0(:knon, :)         coefh = max(coefh, coefh0)
        ycoefm = max(ycoefm, ycoefm0(:knon, :))  
        ycoefh = max(ycoefh, ycoefh0(:knon, :))  
65      END IF      END IF
66    
67      IF (iflag_pbl >= 6) THEN      IF (iflag_pbl >= 6) THEN
68         ! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et         ! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et
69         ! Fr\'ed\'eric Hourdin         ! Fr\'ed\'eric Hourdin
70         yzlay(:knon, 1) = rd * yt(:knon, 1) / (0.5 * (ypaprs(:knon, 1) &         zlay(:, 1) = rd * t(:, 1) / (0.5 * (paprs(:, 1) &
71              + ypplay(:knon, 1))) &              + pplay(:, 1))) * (paprs(:, 1) - pplay(:, 1)) / rg
             * (ypaprs(:knon, 1) - ypplay(:knon, 1)) / rg  
72    
73         DO k = 2, klev         DO k = 2, klev
74            yzlay(:knon, k) = yzlay(:knon, k-1) &            zlay(:, k) = zlay(:, k-1) + rd * 0.5 * (t(:, k-1) + t(:, k)) &
75                 + rd * 0.5 * (yt(1:knon, k-1) + yt(1:knon, k)) &                 / paprs(:, k) * (pplay(:, k-1) - pplay(:, k)) / rg
                / ypaprs(1:knon, k) &  
                * (ypplay(1:knon, k-1) - ypplay(1:knon, k)) / rg  
        END DO  
   
        DO k = 1, klev  
           yteta(1:knon, k) = yt(1:knon, k) * (ypaprs(1:knon, 1) &  
                / ypplay(1:knon, k))**rkappa * (1. + 0.61 * yq(1:knon, k))  
76         END DO         END DO
77    
78         zlev(:knon, 1) = 0.         forall (k = 1:klev) teta(:, k) = t(:, k) &
79         zlev(:knon, klev + 1) = 2. * yzlay(:knon, klev) &              * (paprs(:, 1) / pplay(:, k))**rkappa * (1. + 0.61 * q(:, k))
             - yzlay(:knon, klev - 1)  
80    
81         DO k = 2, klev         zlev(:, 1) = 0.
82            zlev(:knon, k) = 0.5 * (yzlay(:knon, k) + yzlay(:knon, k-1))         zlev(:, klev + 1) = 2. * zlay(:, klev) - zlay(:, klev - 1)
83         END DO         forall (k = 2:klev) zlev(:, k) = 0.5 * (zlay(:, k) + zlay(:, k-1))
84    
85         ustar(:knon) = ustarhb(yu(:knon, 1), yv(:knon, 1), ycdragm(:knon))         CALL yamada4(dtime, zlev, zlay, u, v, teta, q2, coefm, coefh, &
86         CALL yamada4(dtime, rg, zlev(:knon, :), yzlay(:knon, :), &              ustarhb(u(:, 1), v(:, 1), cdragm))
             yu(:knon, :), yv(:knon, :), yteta(:knon, :), yq2(:knon, :), &  
             ycoefm, ycoefh, ustar(:knon))  
87      END IF      END IF
88    
89    end subroutine coef_diff_turb    end subroutine coef_diff_turb

Legend:
Removed from v.250  
changed lines
  Added in v.288

  ViewVC Help
Powered by ViewVC 1.1.21