/[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/phylmd/coef_diff_turb.f revision 254 by guez, Mon Feb 5 10:39:38 2018 UTC trunk/phylmd/Interface_surf/coef_diff_turb.f revision 298 by guez, Thu Jul 26 16:45:51 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, paprs, pplay, u, v, q, t, ts, &    subroutine coef_diff_turb(nsrf, ni, paprs, pplay, u, v, q, t, ts, cdragm, &
8         cdragm, zgeop, coefm, coefh, q2)         zgeop, coefm, coefh, q2)
9    
10      ! Computes coefficients for turbulent diffusion in the atmosphere.      ! 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 comconst, only: dtphys
16      use coefkz_m, only: coefkz      use coefkz_m, only: coefkz
17      use coefkzmin_m, only: coefkzmin      use coefkzmin_m, only: coefkzmin
18      use coefkz2_m, only: coefkz2      use coefkz2_m, only: coefkz2
19      USE conf_phys_m, ONLY: iflag_pbl      USE conf_phys_m, ONLY: iflag_pbl
20      USE dimphy, ONLY: klev      USE dimphy, ONLY: klev
     use nr_util, only: assert  
21      USE suphec_m, ONLY: rd, rg, rkappa      USE suphec_m, ONLY: rd, rg, rkappa
22      use ustarhb_m, only: ustarhb      use ustarhb_m, only: ustarhb
23      use yamada4_m, only: yamada4      use yamada4_m, only: yamada4
24    
     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):: paprs(:, :) ! (knon, klev + 1)      REAL, INTENT(IN):: paprs(:, :) ! (knon, klev + 1)
# Line 39  contains Line 40  contains
40      REAL coefm0(size(ni), 2:klev), coefh0(size(ni), 2:klev) ! (knon, 2:klev)      REAL coefm0(size(ni), 2:klev), coefh0(size(ni), 2:klev) ! (knon, 2:klev)
41      REAL zlay(size(ni), klev), teta(size(ni), klev) ! (knon, klev)      REAL zlay(size(ni), klev), teta(size(ni), klev) ! (knon, klev)
42      real zlev(size(ni), klev + 1) ! (knon, klev + 1)      real zlev(size(ni), klev + 1) ! (knon, klev + 1)
     REAL ustar(size(ni)) ! (knon)  
43      integer k      integer k
44    
45      !-------------------------------------------------------------------------      !-------------------------------------------------------------------------
# Line 48  contains Line 48  contains
48           size(v, 1), size(q, 1), size(t, 1), size(ts), size(cdragm), &           size(v, 1), size(q, 1), size(t, 1), size(ts), size(cdragm), &
49           size(zgeop, 1), size(coefm, 1), size(coefh, 1), size(q2, 1)], &           size(zgeop, 1), size(coefm, 1), size(coefh, 1), size(q2, 1)], &
50           "coef_diff_turb knon")           "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  
51    
52      IF (iflag_pbl >= 6) THEN      IF (iflag_pbl >= 6) THEN
53         ! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et         ! Mellor et Yamada adapt\'e \`a Mars, Richard Fournier et
# Line 76  contains Line 60  contains
60                 / paprs(:, k) * (pplay(:, k-1) - pplay(:, k)) / rg                 / paprs(:, k) * (pplay(:, k-1) - pplay(:, k)) / rg
61         END DO         END DO
62    
63         DO k = 1, klev         forall (k = 1:klev) teta(:, k) = t(:, k) &
64            teta(:, k) = t(:, k) * (paprs(:, 1) / pplay(:, k))**rkappa &              * (paprs(:, 1) / pplay(:, k))**rkappa * (1. + 0.61 * q(:, k))
                * (1. + 0.61 * q(:, k))  
        END DO  
65    
66         zlev(:, 1) = 0.         zlev(:, 1) = 0.
67         zlev(:, klev + 1) = 2. * zlay(:, klev) - zlay(:, klev - 1)         zlev(:, klev + 1) = 2. * zlay(:, klev) - zlay(:, klev - 1)
68           forall (k = 2:klev) zlev(:, k) = 0.5 * (zlay(:, k) + zlay(:, k-1))
69    
70         DO k = 2, klev         CALL yamada4(dtphys, zlev, zlay, u, v, teta, q2, coefm, coefh, &
71            zlev(:, k) = 0.5 * (zlay(:, k) + zlay(:, k-1))              ustarhb(u(:, 1), v(:, 1), cdragm))
72         END DO      else
73           CALL coefkz(nsrf, paprs, pplay, ts, u, v, t, q, zgeop, coefm, coefh)
74         ustar = ustarhb(u(:, 1), v(:, 1), cdragm)  
75         CALL yamada4(dtime, zlev, zlay, u, v, teta, q2, coefm, coefh, ustar)         IF (iflag_pbl == 1) THEN
76              CALL coefkz2(nsrf, paprs, pplay, t, coefm0, coefh0)
77              coefm = max(coefm, coefm0)
78              coefh = max(coefh, coefh0)
79           END IF
80    
81           IF (ok_kzmin) THEN
82              ! Calcul d'une diffusion minimale pour les conditions tres stables
83              CALL coefkzmin(paprs, pplay, u, v, t, q, cdragm, coefh0)
84              coefm = max(coefm, coefh0)
85              coefh = max(coefh, coefh0)
86           END IF
87      END IF      END IF
88    
89    end subroutine coef_diff_turb    end subroutine coef_diff_turb

Legend:
Removed from v.254  
changed lines
  Added in v.298

  ViewVC Help
Powered by ViewVC 1.1.21