/[lmdze]/trunk/Sources/phylmd/CV30_routines/cv30_yield.f
ViewVC logotype

Diff of /trunk/Sources/phylmd/CV30_routines/cv30_yield.f

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

revision 195 by guez, Wed May 18 17:56:44 2016 UTC revision 205 by guez, Tue Jun 21 15:16:03 2016 UTC
# Line 5  module cv30_yield_m Line 5  module cv30_yield_m
5  contains  contains
6    
7    SUBROUTINE cv30_yield(icb, inb, delt, t, rr, u, v, gz, p, ph, h, hp, lv, &    SUBROUTINE cv30_yield(icb, inb, delt, t, rr, u, v, gz, p, ph, h, hp, lv, &
8         cpn, th, ep, clw, m, tp, mp, rp, up, vp, wt, water, evap, b, ment, &         cpn, th, ep, clw, m, tp, mp, qp, up, vp, wt, water, evap, b, ment, &
9         qent, uent, vent, nent, elij, sig, tv, tvp, iflag, precip, VPrecip, &         qent, uent, vent, nent, elij, sig, tv, tvp, iflag, precip, VPrecip, &
10         ft, fr, fu, fv, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)         ft, fr, fu, fv, upwd, dnwd, ma, mike, tls, tps, qcondc)
11    
12      ! Tendencies, precipitation, variables of interface with other      ! Tendencies, precipitation, variables of interface with other
13      ! processes, etc.      ! processes, etc.
14    
15      use conema3_m, only: iflag_clw      use conema3_m, only: iflag_clw
16      use cv30_param_m, only: minorig, nl, sigd      use cv30_param_m, only: minorig, nl, sigd
17      use cv_thermo_m, only: cl, cpd, cpv, grav, rowl, rrd, rrv      use cv_thermo_m, only: rowl
18      USE dimphy, ONLY: klev, klon      USE dimphy, ONLY: klev, klon
19        use SUPHEC_M, only: rg, rcpd, rcw, rcpv, rd, rv
20    
21      ! inputs:      ! inputs:
22      integer, intent(in):: icb(:), inb(:) ! (ncum)  
23        integer, intent(in):: icb(:)
24    
25        integer, intent(in):: inb(:) ! (ncum)
26        ! first model level above the level of neutral buoyancy of the
27        ! parcel (1 <= inb <= nl - 1)
28    
29      real, intent(in):: delt      real, intent(in):: delt
30      real t(klon, klev), rr(klon, klev), u(klon, klev), v(klon, klev)      real, intent(in):: t(klon, klev), rr(klon, klev)
31        real, intent(in):: u(klon, klev), v(klon, klev)
32      real gz(klon, klev)      real gz(klon, klev)
33      real p(klon, klev)      real p(klon, klev)
34      real ph(klon, klev + 1), h(klon, klev), hp(klon, klev)      real ph(klon, klev + 1), h(klon, klev), hp(klon, klev)
35      real lv(klon, klev), cpn(klon, klev)      real, intent(in):: lv(:, :) ! (klon, klev)
36      real th(klon, klev)  
37        real, intent(in):: cpn(:, :) ! (ncum, nl)
38        ! specific heat capacity at constant pressure of humid air, in J K-1 kg-1
39    
40        real, intent(in):: th(:, :) ! (ncum, nl)
41      real ep(klon, klev), clw(klon, klev)      real ep(klon, klev), clw(klon, klev)
42      real m(klon, klev)      real m(klon, klev)
43      real tp(klon, klev)      real tp(klon, klev)
44      real mp(klon, klev), rp(klon, klev), up(klon, klev)  
45        real, intent(in):: mp(:, :) ! (ncum, nl) Mass flux of the
46        ! unsaturated downdraft, defined positive downward, in kg m-2
47        ! s-1. M_p in Emanuel (1991 928).
48    
49        real, intent(in):: qp(:, :), up(:, :) ! (klon, klev)
50      real, intent(in):: vp(:, 2:) ! (ncum, 2:nl)      real, intent(in):: vp(:, 2:) ! (ncum, 2:nl)
51      real, intent(in):: wt(:, :) ! (ncum, nl - 1)      real, intent(in):: wt(:, :) ! (ncum, nl - 1)
52      real, intent(in):: water(:, :), evap(:, :) ! (ncum, nl)      real, intent(in):: water(:, :), evap(:, :) ! (ncum, nl)
# Line 41  contains Line 58  contains
58      real sig(klon, klev)      real sig(klon, klev)
59      real tv(klon, klev), tvp(klon, klev)      real tv(klon, klev), tvp(klon, klev)
60    
     ! input / output:  
     integer iflag(klon)  
   
61      ! outputs:      ! outputs:
62        integer, intent(out):: iflag(:) ! (ncum)
63      real precip(klon)      real precip(klon)
64      real VPrecip(klon, klev + 1)      real VPrecip(klon, klev + 1)
65      real ft(klon, klev), fr(klon, klev), fu(klon, klev), fv(klon, klev)      real ft(klon, klev), fr(klon, klev), fu(klon, klev), fv(klon, klev)
66      real upwd(klon, klev), dnwd(klon, klev)      real upwd(klon, klev), dnwd(klon, klev)
     real dnwd0(klon, klev)  
67      real ma(klon, klev)      real ma(klon, klev)
68      real mike(klon, klev)      real mike(klon, klev)
69      real tls(klon, klev), tps(klon, klev)      real tls(klon, klev), tps(klon, klev)
# Line 58  contains Line 72  contains
72      ! Local:      ! Local:
73      real, parameter:: delta = 0.01 ! interface cloud parameterization      real, parameter:: delta = 0.01 ! interface cloud parameterization
74      integer ncum      integer ncum
75      integer i, k, il, n, j, num1      integer i, k, il, n, j
76      real rat, awat, delti      real awat, delti
77      real ax, bx, cx, dx      real ax, bx, cx, dx
78      real cpinv, rdcp, dpinv      real cpinv, rdcp, dpinv
79      real lvcp(klon, klev)      real lvcp(klon, klev)
# Line 72  contains Line 86  contains
86      !-------------------------------------------------------------      !-------------------------------------------------------------
87    
88      ncum = size(icb)      ncum = size(icb)
89        iflag = 0
90    
91      ! initialization:      ! initialization:
92    
# Line 101  contains Line 116  contains
116         enddo         enddo
117      enddo      enddo
118    
119      ! calculate surface precipitation in mm / day      ! calculate surface precipitation in mm / day
120    
121      do il = 1, ncum      do il = 1, ncum
122         if (ep(il, inb(il)) >= 1e-4) precip(il) = wt(il, 1) * sigd &         if (ep(il, inb(il)) >= 1e-4) precip(il) = wt(il, 1) * sigd &
123              * water(il, 1) * 86400. * 1000. / (rowl * grav)              * water(il, 1) * 86400. * 1000. / (rowl * rg)
124      enddo      enddo
125    
126      ! CALCULATE VERTICAL PROFILE OF PRECIPITATIONs IN kg / m2 / s ===      ! CALCULATE VERTICAL PROFILE OF PRECIPITATIONs IN kg / m2 / s ===
# Line 114  contains Line 129  contains
129      do k = 1, nl - 1      do k = 1, nl - 1
130         do il = 1, ncum         do il = 1, ncum
131            if (k <= inb(il)) VPrecip(il, k) = wt(il, k) * sigd * water(il, k) &            if (k <= inb(il)) VPrecip(il, k) = wt(il, k) * sigd * water(il, k) &
132                 / grav                 / rg
133         end do         end do
134      end do      end do
135    
136      ! calculate tendencies of lowest level potential temperature      ! calculate tendencies of lowest level potential temperature
137      ! and mixing ratio      ! and mixing ratio
138    
139      do il = 1, ncum      do il = 1, ncum
140         work(il) = 1.0 / (ph(il, 1) - ph(il, 2))         work(il) = 1.0 / (ph(il, 1) - ph(il, 2))
# Line 133  contains Line 148  contains
148      enddo      enddo
149    
150      do il = 1, ncum      do il = 1, ncum
151         ! Consist vect:         if (0.01 * rg * work(il) * am(il) >= delti) iflag(il) = 1
        if (0.01 * grav * work(il) * am(il) >= delti) iflag(il) = 1  
   
        ft(il, 1) = 0.01 * grav * work(il) * am(il) * (t(il, 2) - t(il, 1) &  
             + (gz(il, 2) - gz(il, 1)) / cpn(il, 1))  
   
        ft(il, 1) = ft(il, 1) - 0.5 * lvcp(il, 1) * sigd * (evap(il, 1) &  
             + evap(il, 2))  
   
        ft(il, 1) = ft(il, 1) - 0.009 * grav * sigd * mp(il, 2) &  
             * t(il, 1) * b(il, 1) * work(il)  
152    
153         ft(il, 1) = ft(il, 1) + 0.01 * sigd * wt(il, 1) * (cl - cpd) &         ft(il, 1) = 0.01 * rg * work(il) * am(il) * (t(il, 2) - t(il, 1) &
154              * water(il, 2) * (t(il, 2) - t(il, 1)) * work(il) / cpn(il, 1)              + (gz(il, 2) - gz(il, 1)) / cpn(il, 1)) - 0.5 * lvcp(il, 1) &
155                * sigd * (evap(il, 1) + evap(il, 2)) - 0.009 * rg * sigd &
156         !jyg1 Correction pour mieux conserver l'eau (conformite avec CONVECT4.3)              * mp(il, 2) * t(il, 1) * b(il, 1) * work(il) + 0.01 * sigd &
157         ! (sb: pour l'instant, on ne fait que le chgt concernant grav, pas evap)              * wt(il, 1) * (rcw - rcpd) * water(il, 2) * (t(il, 2) - t(il, 1)) &
158         fr(il, 1) = 0.01 * grav * mp(il, 2) * (rp(il, 2) - rr(il, 1)) &              * work(il) / cpn(il, 1)
159    
160           ! jyg Correction pour mieux conserver l'eau (conformite avec CONVECT4.3)
161           ! (sb: pour l'instant, on ne fait que le chgt concernant rg, pas evap)
162           fr(il, 1) = 0.01 * rg * mp(il, 2) * (qp(il, 2) - rr(il, 1)) &
163              * work(il) + sigd * 0.5 * (evap(il, 1) + evap(il, 2))              * work(il) + sigd * 0.5 * (evap(il, 1) + evap(il, 2))
164         ! + tard : + sigd * evap(il, 1)         ! + tard : + sigd * evap(il, 1)
165    
166         fr(il, 1) = fr(il, 1) + 0.01 * grav * am(il) * (rr(il, 2) - rr(il, 1)) &         fr(il, 1) = fr(il, 1) + 0.01 * rg * am(il) * (rr(il, 2) - rr(il, 1)) &
167              * work(il)              * work(il)
168    
169         fu(il, 1) = fu(il, 1) + 0.01 * grav * work(il) * (mp(il, 2) &         fu(il, 1) = fu(il, 1) + 0.01 * rg * work(il) * (mp(il, 2) &
170              * (up(il, 2) - u(il, 1)) + am(il) * (u(il, 2) - u(il, 1)))              * (up(il, 2) - u(il, 1)) + am(il) * (u(il, 2) - u(il, 1)))
171         fv(il, 1) = fv(il, 1) + 0.01 * grav * work(il) * (mp(il, 2) &         fv(il, 1) = fv(il, 1) + 0.01 * rg * work(il) * (mp(il, 2) &
172              * (vp(il, 2) - v(il, 1)) + am(il) * (v(il, 2) - v(il, 1)))              * (vp(il, 2) - v(il, 1)) + am(il) * (v(il, 2) - v(il, 1)))
173      enddo ! il      enddo
174    
175      do j = 2, nl      do j = 2, nl
176         do il = 1, ncum         do il = 1, ncum
177            if (j <= inb(il)) then            if (j <= inb(il)) then
178               fr(il, 1) = fr(il, 1) + 0.01 * grav * work(il) * ment(il, j, 1) &               fr(il, 1) = fr(il, 1) + 0.01 * rg * work(il) * ment(il, j, 1) &
179                    * (qent(il, j, 1) - rr(il, 1))                    * (qent(il, j, 1) - rr(il, 1))
180               fu(il, 1) = fu(il, 1) + 0.01 * grav * work(il) * ment(il, j, 1) &               fu(il, 1) = fu(il, 1) + 0.01 * rg * work(il) * ment(il, j, 1) &
181                    * (uent(il, j, 1) - u(il, 1))                    * (uent(il, j, 1) - u(il, 1))
182               fv(il, 1) = fv(il, 1) + 0.01 * grav * work(il) * ment(il, j, 1) &               fv(il, 1) = fv(il, 1) + 0.01 * rg * work(il) * ment(il, j, 1) &
183                    * (vent(il, j, 1) - v(il, 1))                    * (vent(il, j, 1) - v(il, 1))
184            endif            endif
185         enddo         enddo
186      enddo      enddo
187    
188      ! calculate tendencies of potential temperature and mixing ratio      ! calculate tendencies of potential temperature and mixing ratio
189      ! at levels above the lowest level      ! at levels above the lowest level
190    
191      ! first find the net saturated updraft and downdraft mass fluxes      ! first find the net saturated updraft and downdraft mass fluxes
192      ! through each level      ! through each level
193    
194      loop_i: do i = 2, nl - 1      loop_i: do i = 2, nl - 1
195         num1 = 0         if (any(inb >= i)) then
   
        do il = 1, ncum  
           if (i <= inb(il)) num1 = num1 + 1  
        enddo  
   
        if (num1 > 0) then  
196            amp1(:ncum) = 0.            amp1(:ncum) = 0.
197            ad(:ncum) = 0.            ad(:ncum) = 0.
198    
# Line 226  contains Line 229  contains
229                  dpinv = 1.0 / (ph(il, i) - ph(il, i + 1))                  dpinv = 1.0 / (ph(il, i) - ph(il, i + 1))
230                  cpinv = 1.0 / cpn(il, i)                  cpinv = 1.0 / cpn(il, i)
231    
232                  ! Vecto:                  if (0.01 * rg * dpinv * amp1(il) >= delti) iflag(il) = 1
                 if (0.01 * grav * dpinv * amp1(il) >= delti) iflag(il) = 1  
233    
234                  ft(il, i) = 0.01 * grav * dpinv * (amp1(il) * (t(il, i + 1) &                  ft(il, i) = 0.01 * rg * dpinv * (amp1(il) * (t(il, i + 1) &
235                       - t(il, i) + (gz(il, i + 1) - gz(il, i)) * cpinv) &                       - t(il, i) + (gz(il, i + 1) - gz(il, i)) * cpinv) &
236                       - ad(il) * (t(il, i) - t(il, i - 1) + (gz(il, i) &                       - ad(il) * (t(il, i) - t(il, i - 1) + (gz(il, i) &
237                       - gz(il, i - 1)) * cpinv)) - 0.5 * sigd * lvcp(il, i) &                       - gz(il, i - 1)) * cpinv)) - 0.5 * sigd * lvcp(il, i) &
238                       * (evap(il, i) + evap(il, i + 1))                       * (evap(il, i) + evap(il, i + 1)) - 0.009 * rg * sigd &
239                  rat = cpn(il, i - 1) * cpinv                       * (mp(il, i + 1) * t(il, i) * b(il, i) - mp(il, i) &
240                  ft(il, i) = ft(il, i) - 0.009 * grav * sigd * (mp(il, i + 1) &                       * t(il, i - 1) * cpn(il, i - 1) * cpinv * b(il, i - 1)) &
241                       * t(il, i) * b(il, i) - mp(il, i) * t(il, i - 1) * rat &                       * dpinv + 0.01 * rg * dpinv * ment(il, i, i) &
242                       * b(il, i - 1)) * dpinv                       * (hp(il, i) - h(il, i) + t(il, i) * (rcpv - rcpd) &
243                  ft(il, i) = ft(il, i) + 0.01 * grav * dpinv * ment(il, i, i) &                       * (rr(il, i) - qent(il, i, i))) * cpinv + 0.01 * sigd &
244                       * (hp(il, i) - h(il, i) + t(il, i) * (cpv - cpd) &                       * wt(il, i) * (rcw - rcpd) * water(il, i + 1) &
245                       * (rr(il, i) - qent(il, i, i))) * cpinv                       * (t(il, i + 1) - t(il, i)) * dpinv * cpinv
246                    fr(il, i) = 0.01 * rg * dpinv * (amp1(il) * (rr(il, i + 1) &
                 ft(il, i) = ft(il, i) + 0.01 * sigd * wt(il, i) * (cl - cpd) &  
                      * water(il, i + 1) * (t(il, i + 1) - t(il, i)) * dpinv &  
                      * cpinv  
   
                 fr(il, i) = 0.01 * grav * dpinv * (amp1(il) * (rr(il, i + 1) &  
247                       - rr(il, i)) - ad(il) * (rr(il, i) - rr(il, i - 1)))                       - rr(il, i)) - ad(il) * (rr(il, i) - rr(il, i - 1)))
248                  fu(il, i) = fu(il, i) + 0.01 * grav * dpinv * (amp1(il) &                  fu(il, i) = fu(il, i) + 0.01 * rg * dpinv * (amp1(il) &
249                       * (u(il, i + 1) - u(il, i)) - ad(il) * (u(il, i) &                       * (u(il, i + 1) - u(il, i)) - ad(il) * (u(il, i) &
250                       - u(il, i - 1)))                       - u(il, i - 1)))
251                  fv(il, i) = fv(il, i) + 0.01 * grav * dpinv * (amp1(il) &                  fv(il, i) = fv(il, i) + 0.01 * rg * dpinv * (amp1(il) &
252                       * (v(il, i + 1) - v(il, i)) - ad(il) * (v(il, i) &                       * (v(il, i + 1) - v(il, i)) - ad(il) * (v(il, i) &
253                       - v(il, i - 1)))                       - v(il, i - 1)))
254               endif               endif
# Line 266  contains Line 263  contains
263                     awat = elij(il, k, i) - (1. - ep(il, i)) * clw(il, i)                     awat = elij(il, k, i) - (1. - ep(il, i)) * clw(il, i)
264                     awat = amax1(awat, 0.0)                     awat = amax1(awat, 0.0)
265    
266                     fr(il, i) = fr(il, i) + 0.01 * grav * dpinv &                     fr(il, i) = fr(il, i) + 0.01 * rg * dpinv &
267                          * ment(il, k, i) * (qent(il, k, i) - awat - rr(il, i))                          * ment(il, k, i) * (qent(il, k, i) - awat - rr(il, i))
268                     fu(il, i) = fu(il, i) + 0.01 * grav * dpinv &                     fu(il, i) = fu(il, i) + 0.01 * rg * dpinv &
269                          * ment(il, k, i) * (uent(il, k, i) - u(il, i))                          * ment(il, k, i) * (uent(il, k, i) - u(il, i))
270                     fv(il, i) = fv(il, i) + 0.01 * grav * dpinv &                     fv(il, i) = fv(il, i) + 0.01 * rg * dpinv &
271                          * ment(il, k, i) * (vent(il, k, i) - v(il, i))                          * ment(il, k, i) * (vent(il, k, i) - v(il, i))
272    
273                     ! (saturated updrafts resulting from mixing)                     ! (saturated updrafts resulting from mixing)
# Line 286  contains Line 283  contains
283                     dpinv = 1.0 / (ph(il, i) - ph(il, i + 1))                     dpinv = 1.0 / (ph(il, i) - ph(il, i + 1))
284                     cpinv = 1.0 / cpn(il, i)                     cpinv = 1.0 / cpn(il, i)
285    
286                     fr(il, i) = fr(il, i) + 0.01 * grav * dpinv &                     fr(il, i) = fr(il, i) + 0.01 * rg * dpinv &
287                          * ment(il, k, i) * (qent(il, k, i) - rr(il, i))                          * ment(il, k, i) * (qent(il, k, i) - rr(il, i))
288                     fu(il, i) = fu(il, i) + 0.01 * grav * dpinv &                     fu(il, i) = fu(il, i) + 0.01 * rg * dpinv &
289                          * ment(il, k, i) * (uent(il, k, i) - u(il, i))                          * ment(il, k, i) * (uent(il, k, i) - u(il, i))
290                     fv(il, i) = fv(il, i) + 0.01 * grav * dpinv &                     fv(il, i) = fv(il, i) + 0.01 * rg * dpinv &
291                          * ment(il, k, i) * (vent(il, k, i) - v(il, i))                          * ment(il, k, i) * (vent(il, k, i) - v(il, i))
292                  endif                  endif
293               end do               end do
# Line 304  contains Line 301  contains
301                  ! sb: on ne fait pas encore la correction permettant de mieux                  ! sb: on ne fait pas encore la correction permettant de mieux
302                  ! conserver l'eau:                  ! conserver l'eau:
303                  fr(il, i) = fr(il, i) + 0.5 * sigd * (evap(il, i) &                  fr(il, i) = fr(il, i) + 0.5 * sigd * (evap(il, i) &
304                       + evap(il, i + 1)) + 0.01 * grav * (mp(il, i + 1) &                       + evap(il, i + 1)) + 0.01 * rg * (mp(il, i + 1) &
305                       * (rp(il, i + 1) - rr(il, i)) - mp(il, i) * (rp(il, i) &                       * (qp(il, i + 1) - rr(il, i)) - mp(il, i) * (qp(il, i) &
306                       - rr(il, i - 1))) * dpinv                       - rr(il, i - 1))) * dpinv
307    
308                  fu(il, i) = fu(il, i) + 0.01 * grav * (mp(il, i + 1) &                  fu(il, i) = fu(il, i) + 0.01 * rg * (mp(il, i + 1) &
309                       * (up(il, i + 1) - u(il, i)) - mp(il, i) * (up(il, i) &                       * (up(il, i + 1) - u(il, i)) - mp(il, i) * (up(il, i) &
310                       - u(il, i - 1))) * dpinv                       - u(il, i - 1))) * dpinv
311                  fv(il, i) = fv(il, i) + 0.01 * grav * (mp(il, i + 1) &                  fv(il, i) = fv(il, i) + 0.01 * rg * (mp(il, i + 1) &
312                       * (vp(il, i + 1) - v(il, i)) - mp(il, i) * (vp(il, i) &                       * (vp(il, i + 1) - v(il, i)) - mp(il, i) * (vp(il, i) &
313                       - v(il, i - 1))) * dpinv                       - v(il, i - 1))) * dpinv
314               endif               endif
# Line 345  contains Line 342  contains
342         end if         end if
343      end do loop_i      end do loop_i
344    
345      ! move the detrainment at level inb down to level inb - 1      ! move the detrainment at level inb down to level inb - 1
346      ! in such a way as to preserve the vertically      ! in such a way as to preserve the vertically
347      ! integrated enthalpy and water tendencies      ! integrated enthalpy and water tendencies
348    
349      do il = 1, ncum      do il = 1, ncum
350         ax = 0.1 * ment(il, inb(il), inb(il)) * (hp(il, inb(il)) &         ax = 0.1 * ment(il, inb(il), inb(il)) * (hp(il, inb(il)) &
351              - h(il, inb(il)) + t(il, inb(il)) * (cpv - cpd) &              - h(il, inb(il)) + t(il, inb(il)) * (rcpv - rcpd) &
352              * (rr(il, inb(il)) - qent(il, inb(il), inb(il)))) &              * (rr(il, inb(il)) - qent(il, inb(il), inb(il)))) &
353              / (cpn(il, inb(il)) * (ph(il, inb(il)) - ph(il, inb(il) + 1)))              / (cpn(il, inb(il)) * (ph(il, inb(il)) - ph(il, inb(il) + 1)))
354         ft(il, inb(il)) = ft(il, inb(il)) - ax         ft(il, inb(il)) = ft(il, inb(il)) - ax
# Line 382  contains Line 379  contains
379    
380      end do      end do
381    
382      ! homoginize tendencies below cloud base      ! homoginize tendencies below cloud base
383    
384      do il = 1, ncum      do il = 1, ncum
385         asum(il) = 0.0         asum(il) = 0.0
# Line 395  contains Line 392  contains
392         do il = 1, ncum         do il = 1, ncum
393            if (i <= (icb(il) - 1)) then            if (i <= (icb(il) - 1)) then
394               asum(il) = asum(il) + ft(il, i) * (ph(il, i) - ph(il, i + 1))               asum(il) = asum(il) + ft(il, i) * (ph(il, i) - ph(il, i + 1))
395               bsum(il) = bsum(il) + fr(il, i) * (lv(il, i) + (cl - cpd) &               bsum(il) = bsum(il) + fr(il, i) * (lv(il, i) + (rcw - rcpd) &
396                    * (t(il, i) - t(il, 1))) * (ph(il, i) - ph(il, i + 1))                    * (t(il, i) - t(il, 1))) * (ph(il, i) - ph(il, i + 1))
397               csum(il) = csum(il) + (lv(il, i) + (cl - cpd) * (t(il, i) &               csum(il) = csum(il) + (lv(il, i) + (rcw - rcpd) * (t(il, i) &
398                    - t(il, 1))) * (ph(il, i) - ph(il, i + 1))                    - t(il, 1))) * (ph(il, i) - ph(il, i + 1))
399               dsum(il) = dsum(il) + t(il, i) * (ph(il, i) - ph(il, i + 1)) &               dsum(il) = dsum(il) + t(il, i) * (ph(il, i) - ph(il, i + 1)) &
400                    / th(il, i)                    / th(il, i)
# Line 414  contains Line 411  contains
411         enddo         enddo
412      enddo      enddo
413    
414      ! reset counter and return      ! reset counter and return
415    
416      do il = 1, ncum      do il = 1, ncum
417         sig(il, klev) = 2.0         sig(il, klev) = 2.0
# Line 429  contains Line 426  contains
426    
427      do i = 1, nl      do i = 1, nl
428         do il = 1, ncum         do il = 1, ncum
           dnwd0(il, i) = - mp(il, i)  
        enddo  
     enddo  
     do i = nl + 1, klev  
        do il = 1, ncum  
           dnwd0(il, i) = 0.  
        enddo  
     enddo  
   
     do i = 1, nl  
        do il = 1, ncum  
429            if (i >= icb(il) .and. i <= inb(il)) then            if (i >= icb(il) .and. i <= inb(il)) then
430               upwd(il, i) = 0.0               upwd(il, i) = 0.0
431               dnwd(il, i) = 0.0               dnwd(il, i) = 0.0
# Line 480  contains Line 466  contains
466         enddo         enddo
467      enddo      enddo
468    
469      !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc      ! D\'etermination de la variation de flux ascendant entre
470      ! determination de la variation de flux ascendant entre      ! deux niveaux non dilu\'es mike
     ! deux niveau non dilue mike  
     !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc  
471    
472      do i = 1, nl      do i = 1, nl
473         do il = 1, ncum         do il = 1, ncum
# Line 525  contains Line 509  contains
509         enddo         enddo
510      enddo      enddo
511    
512      !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc      ! icb repr\'esente le niveau o\`u se trouve la base du nuage, et
513      ! icb represente de niveau ou se trouve la      ! inb le sommet du nuage
     ! base du nuage, et inb le top du nuage  
     !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc  
514    
515      do i = 1, klev      do i = 1, klev
516         DO il = 1, ncum         DO il = 1, ncum
517            rdcp = (rrd * (1. - rr(il, i)) - rr(il, i) * rrv) &            rdcp = (rd * (1. - rr(il, i)) - rr(il, i) * rv) &
518                 / (cpd * (1. - rr(il, i)) + rr(il, i) * cpv)                 / (rcpd * (1. - rr(il, i)) + rr(il, i) * rcpv)
519            tls(il, i) = t(il, i) * (1000.0 / p(il, i))**rdcp            tls(il, i) = t(il, i) * (1000.0 / p(il, i))**rdcp
520            tps(il, i) = tp(il, i)            tps(il, i) = tp(il, i)
521         end DO         end DO
522      enddo      enddo
523    
524      ! diagnose the in-cloud mixing ratio      ! Diagnose the in-cloud mixing ratio of condensed water
     ! of condensed water  
     !  
525    
526      do i = 1, klev      do i = 1, klev
527         do il = 1, ncum         do il = 1, ncum
# Line 567  contains Line 547  contains
547            do il = 1, ncum            do il = 1, ncum
548               if (i >= icb(il) .and. i <= (inb(il) - 1) &               if (i >= icb(il) .and. i <= (inb(il) - 1) &
549                    .and. j >= icb(il)) then                    .and. j >= icb(il)) then
550                  sax(il, i) = sax(il, i) + rrd * (tvp(il, j) - tv(il, j)) &                  sax(il, i) = sax(il, i) + rd * (tvp(il, j) - tv(il, j)) &
551                       * (ph(il, j) - ph(il, j + 1)) / p(il, j)                       * (ph(il, j) - ph(il, j + 1)) / p(il, j)
552               endif               endif
553            enddo            enddo
# Line 585  contains Line 565  contains
565    
566      do i = 1, nl      do i = 1, nl
567         do il = 1, ncum         do il = 1, ncum
568            if (wa(il, i) > 0.0) siga(il, i) = mac(il, i) / wa(il, i) * rrd &            if (wa(il, i) > 0.0) siga(il, i) = mac(il, i) / wa(il, i) * rd &
569                 * tvp(il, i) / p(il, i) / 100. / delta                 * tvp(il, i) / p(il, i) / 100. / delta
570            siga(il, i) = min(siga(il, i), 1.0)            siga(il, i) = min(siga(il, i), 1.0)
571    

Legend:
Removed from v.195  
changed lines
  Added in v.205

  ViewVC Help
Powered by ViewVC 1.1.21