/[lmdze]/trunk/Sources/phylmd/cv_driver.f
ViewVC logotype

Diff of /trunk/Sources/phylmd/cv_driver.f

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

revision 190 by guez, Thu Apr 14 15:15:56 2016 UTC revision 201 by guez, Mon Jun 6 17:42:15 2016 UTC
# Line 5  module cv_driver_m Line 5  module cv_driver_m
5  contains  contains
6    
7    SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, fq1, fu1, &    SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, fq1, fu1, &
8         fv1, precip1, VPrecip1, sig1, w01, icb1, inb1, delt, Ma1, upwd1, dnwd1, &         fv1, precip1, VPrecip1, sig1, w01, icb1, inb1, Ma1, upwd1, dnwd1, &
9         dnwd01, qcondc1, cape1, da1, phi1, mp1)         dnwd01, qcondc1, cape1, da1, phi1, mp1)
10    
11      ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17      ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
# Line 14  contains Line 14  contains
14    
15      ! Several modules corresponding to different physical processes      ! Several modules corresponding to different physical processes
16    
17        use comconst, only: dtphys
18      use cv30_closure_m, only: cv30_closure      use cv30_closure_m, only: cv30_closure
19      use cv30_compress_m, only: cv30_compress      use cv30_compress_m, only: cv30_compress
20      use cv30_feed_m, only: cv30_feed      use cv30_feed_m, only: cv30_feed
# Line 23  contains Line 24  contains
24      use cv30_tracer_m, only: cv30_tracer      use cv30_tracer_m, only: cv30_tracer
25      use cv30_trigger_m, only: cv30_trigger      use cv30_trigger_m, only: cv30_trigger
26      use cv30_uncompress_m, only: cv30_uncompress      use cv30_uncompress_m, only: cv30_uncompress
27        use cv30_undilute1_m, only: cv30_undilute1
28      use cv30_undilute2_m, only: cv30_undilute2      use cv30_undilute2_m, only: cv30_undilute2
29      use cv30_unsat_m, only: cv30_unsat      use cv30_unsat_m, only: cv30_unsat
30      use cv30_yield_m, only: cv30_yield      use cv30_yield_m, only: cv30_yield
     use cv_thermo_m, only: cv_thermo  
31      USE dimphy, ONLY: klev, klon      USE dimphy, ONLY: klev, klon
32    
33      real, intent(in):: t1(klon, klev) ! temperature (K)      real, intent(in):: t1(klon, klev) ! temperature, in K
34      real, intent(in):: q1(klon, klev) ! specific humidity      real, intent(in):: q1(klon, klev) ! specific humidity
35      real, intent(in):: qs1(klon, klev) ! saturation specific humidity      real, intent(in):: qs1(klon, klev) ! saturation specific humidity
36    
37      real, intent(in):: u1(klon, klev), v1(klon, klev)      real, intent(in):: u1(klon, klev), v1(klon, klev)
38      ! zonal wind and meridional velocity (m/s)      ! zonal wind and meridional velocity (m/s)
39    
40      real, intent(in):: p1(klon, klev) ! full level pressure (hPa)      real, intent(in):: p1(klon, klev) ! full level pressure, in hPa
41    
42      real, intent(in):: ph1(klon, klev + 1)      real, intent(in):: ph1(klon, klev + 1)
43      ! Half level pressure (hPa). These pressures are defined at levels      ! Half level pressure, in hPa. These pressures are defined at levels
44      ! intermediate between those of P1, T1, Q1 and QS1. The first      ! intermediate between those of P1, T1, Q1 and QS1. The first
45      ! value of PH should be greater than (i.e. at a lower level than)      ! value of PH should be greater than (i.e. at a lower level than)
46      ! the first value of the array P1.      ! the first value of the array P1.
47    
48      integer, intent(out):: iflag1(klon)      integer, intent(out):: iflag1(:) ! (klon)
49      ! Flag for Emanuel conditions.      ! Flag for Emanuel conditions.
50    
51      ! 0: Moist convection occurs.      ! 0: Moist convection occurs.
# Line 57  contains Line 58  contains
58    
59      ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.      ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
60    
61      ! 4: No moist convection; atmosphere is not unstable      ! 4: No moist convection; atmosphere is not unstable.
62    
63      ! 6: No moist convection because ihmin le minorig.      ! 6: No moist convection because ihmin <= minorig.
64    
65      ! 7: No moist convection because unreasonable parcel level      ! 7: No moist convection because unreasonable parcel level
66      ! temperature or specific humidity.      ! temperature or specific humidity.
67    
68      ! 8: No moist convection: lifted condensation level is above the      ! 8: No moist convection: lifted condensation level is above the
69      ! 200 mb level.      ! 200 mbar level.
70    
71      ! 9: No moist convection: cloud base is higher then the level NL-1.      ! 9: No moist convection: cloud base is higher than the level NL-1.
72    
73      real, intent(out):: ft1(klon, klev) ! temperature tendency (K/s)      real, intent(out):: ft1(klon, klev) ! temperature tendency (K/s)
74      real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1)      real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1)
# Line 87  contains Line 88  contains
88    
89      integer, intent(out):: icb1(klon)      integer, intent(out):: icb1(klon)
90      integer, intent(inout):: inb1(klon)      integer, intent(inout):: inb1(klon)
     real, intent(in):: delt ! the model time step (sec) between calls  
   
91      real, intent(out):: Ma1(klon, klev) ! mass flux of adiabatic updraft      real, intent(out):: Ma1(klon, klev) ! mass flux of adiabatic updraft
92    
93      real, intent(out):: upwd1(klon, klev)      real, intent(out):: upwd1(klon, klev)
# Line 106  contains Line 105  contains
105    
106      ! Local:      ! Local:
107    
108      real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)      real da(klon, klev), phi(klon, klev, klev)
109        real, allocatable:: mp(:, :) ! (ncum, nl)
110      integer i, k, il      integer i, k, il
     integer icbmax  
     integer nk1(klon)  
111      integer icbs1(klon)      integer icbs1(klon)
112      real plcl1(klon)      real plcl1(klon)
113      real tnk1(klon)      real tnk1(klon)
# Line 117  contains Line 115  contains
115      real gznk1(klon)      real gznk1(klon)
116      real pbase1(klon)      real pbase1(klon)
117      real buoybase1(klon)      real buoybase1(klon)
118      real lv1(klon, klev)  
119      real cpn1(klon, klev)      real lv1(klon, nl)
120        ! specific latent heat of vaporization of water, in J kg-1
121    
122        real cpn1(klon, nl)
123        ! specific heat capacity at constant pressure of humid air, in J K-1 kg-1
124    
125      real tv1(klon, klev)      real tv1(klon, klev)
126      real gz1(klon, klev)      real gz1(klon, klev)
127      real hm1(klon, klev)      real hm1(klon, klev)
# Line 126  contains Line 129  contains
129      real tp1(klon, klev)      real tp1(klon, klev)
130      real tvp1(klon, klev)      real tvp1(klon, klev)
131      real clw1(klon, klev)      real clw1(klon, klev)
132      real th1(klon, klev)      real th1(klon, nl) ! potential temperature, in K
133      integer ncum      integer ncum
134    
135      ! Compressed fields:      ! Compressed fields:
136      integer idcum(klon)      integer, allocatable:: idcum(:), iflag(:) ! (ncum)
137      integer iflag(klon), nk(klon), icb(klon)      integer, allocatable:: icb(:) ! (ncum)
138      integer nent(klon, klev)      integer nent(klon, klev)
139      integer icbs(klon)      integer icbs(klon)
140      integer inb(klon)  
141      real plcl(klon), tnk(klon), qnk(klon), gznk(klon)      integer, allocatable:: inb(:) ! (ncum)
142        ! first model level above the level of neutral buoyancy of the
143        ! parcel (1 <= inb <= nl - 1)
144    
145        real, allocatable:: plcl(:) ! (ncum)
146        real tnk(klon), qnk(klon), gznk(klon)
147      real t(klon, klev), q(klon, klev), qs(klon, klev)      real t(klon, klev), q(klon, klev), qs(klon, klev)
148      real u(klon, klev), v(klon, klev)      real u(klon, klev), v(klon, klev)
149      real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)      real gz(klon, klev), h(klon, klev)
150      real p(klon, klev), ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)  
151        real, allocatable:: lv(:, :) ! (ncum, nl)
152        ! specific latent heat of vaporization of water, in J kg-1
153    
154        real, allocatable:: cpn(:, :) ! (ncum, nl)
155        ! specific heat capacity at constant pressure of humid air, in J K-1 kg-1
156    
157        real p(klon, klev) ! pressure at full level, in hPa
158        real ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)
159      real clw(klon, klev)      real clw(klon, klev)
160      real pbase(klon), buoybase(klon), th(klon, klev)      real pbase(klon), buoybase(klon)
161        real, allocatable:: th(:, :) ! (ncum, nl)
162      real tvp(klon, klev)      real tvp(klon, klev)
163      real sig(klon, klev), w0(klon, klev)      real sig(klon, klev), w0(klon, klev)
164      real hp(klon, klev), ep(klon, klev), sigp(klon, klev)      real hp(klon, klev), ep(klon, klev)
165      real buoy(klon, klev)      real buoy(klon, klev)
166      real cape(klon)      real cape(klon)
167      real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)      real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
# Line 152  contains Line 169  contains
169      real ments(klon, klev, klev), qents(klon, klev, klev)      real ments(klon, klev, klev), qents(klon, klev, klev)
170      real sij(klon, klev, klev), elij(klon, klev, klev)      real sij(klon, klev, klev), elij(klon, klev, klev)
171      real qp(klon, klev), up(klon, klev), vp(klon, klev)      real qp(klon, klev), up(klon, klev), vp(klon, klev)
172      real wt(klon, klev), water(klon, klev), evap(klon, klev)      real wt(klon, klev), water(klon, klev)
173        real, allocatable:: evap(:, :) ! (ncum, nl)
174      real, allocatable:: b(:, :) ! (ncum, nl - 1)      real, allocatable:: b(:, :) ! (ncum, nl - 1)
175      real ft(klon, klev), fq(klon, klev)      real ft(klon, klev), fq(klon, klev)
176      real fu(klon, klev), fv(klon, klev)      real fu(klon, klev), fv(klon, klev)
# Line 166  contains Line 184  contains
184      !-------------------------------------------------------------------      !-------------------------------------------------------------------
185    
186      ! SET CONSTANTS AND PARAMETERS      ! SET CONSTANTS AND PARAMETERS
187        CALL cv30_param
     ! set thermodynamical constants:  
     CALL cv_thermo  
   
     ! set convect parameters  
     ! includes microphysical parameters and parameters that  
     ! control the rate of approach to quasi-equilibrium)  
     ! (common cvparam)  
     CALL cv30_param(delt)  
188    
189      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS      ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
190    
# Line 198  contains Line 208  contains
208         end do         end do
209      end do      end do
210    
211      do i = 1, klon      precip1 = 0.
212         precip1(i) = 0.      cape1 = 0.
213         iflag1(i) = 0      VPrecip1(:, klev + 1) = 0.
        cape1(i) = 0.  
        VPrecip1(i, klev + 1) = 0.  
     end do  
214    
215      do il = 1, klon      do il = 1, klon
216         sig1(il, klev) = sig1(il, klev) + 1.         sig1(il, klev) = sig1(il, klev) + 1.
217         sig1(il, klev) = min(sig1(il, klev), 12.1)         sig1(il, klev) = min(sig1(il, klev), 12.1)
218      enddo      enddo
219    
220      ! CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY & STATIC ENERGY      CALL cv30_prelim(t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, hm1, th1)
221      CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &      CALL cv30_feed(t1, q1, qs1, p1, ph1, gz1, icb1, iflag1, tnk1, qnk1, &
222           gz1, h1, hm1, th1)           gznk1, plcl1)
223        CALL cv30_undilute1(t1, q1, qs1, gz1, plcl1, p1, icb1, tp1, tvp1, clw1, &
224      ! CONVECTIVE FEED           icbs1)
225      CALL cv30_feed(klon, klev, t1, q1, qs1, p1, ph1, gz1, nk1, icb1, &      CALL cv30_trigger(icb1, plcl1, p1, th1, tv1, tvp1, pbase1, buoybase1, &
226           icbmax, iflag1, tnk1, qnk1, gznk1, plcl1) ! klev->na           iflag1, sig1, w01)
227    
228      CALL cv30_undilute1(klon, klev, t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, &      ncum = count(iflag1 == 0)
          tp1, tvp1, clw1, icbs1) ! klev->na  
   
     ! TRIGGERING  
     CALL cv30_trigger(klon, klev, icb1, plcl1, p1, th1, tv1, tvp1, pbase1, &  
          buoybase1, iflag1, sig1, w01) ! klev->na  
   
     ! Moist convective adjustment is necessary  
   
     ncum = 0  
     do i = 1, klon  
        if (iflag1(i) == 0) then  
           ncum = ncum + 1  
           idcum(ncum) = i  
        endif  
     end do  
229    
230      IF (ncum > 0) THEN      IF (ncum > 0) THEN
231         allocate(b(ncum, nl - 1))         ! Moist convective adjustment is necessary
232         CALL cv30_compress(ncum, iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, &         allocate(idcum(ncum), plcl(ncum), inb(ncum))
233           allocate(b(ncum, nl - 1), evap(ncum, nl), icb(ncum), iflag(ncum))
234           allocate(th(ncum, nl), lv(ncum, nl), cpn(ncum, nl), mp(ncum, nl))
235           idcum = pack((/(i, i = 1, klon)/), iflag1 == 0)
236           CALL cv30_compress(idcum, iflag1, icb1, icbs1, plcl1, tnk1, qnk1, &
237              gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, &              gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, &
238              cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, nk, icb, &              cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, icb, icbs, plcl, &
239              icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, &              tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, &
240              th, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, w0)              cpn, p, ph, tv, tp, tvp, clw, sig, w0)
241         CALL cv30_undilute2(ncum, icb, icbs, nk, tnk, qnk, gznk, t, qs, gz, p, &         CALL cv30_undilute2(icb, icbs(:ncum), tnk, qnk, gznk, t, qs, gz, p, h, &
242              h, tv, lv, pbase, buoybase, plcl, inb(:ncum), tp, tvp, clw, hp, &              tv, lv, pbase(:ncum), buoybase(:ncum), plcl, inb, tp, tvp, &
243              ep, sigp, buoy)              clw, hp, ep, buoy)
244           CALL cv30_closure(icb, inb, pbase, p, ph(:ncum, :), tv, buoy, &
245         ! CLOSURE              sig, w0, cape, m)
246         CALL cv30_closure(klon, ncum, klev, icb, inb, pbase, p, ph, tv, &         CALL cv30_mixing(icb, inb, t, q, qs, u, v, h, lv, &
247              buoy, sig, w0, cape, m) ! na->klev              hp, ep, clw, m, sig, ment, qent, uent, vent, nent, sij, elij, &
248                ments, qents)
249         ! MIXING         CALL cv30_unsat(icb, inb, t(:ncum, :nl), q(:ncum, :nl), &
250         CALL cv30_mixing(klon, ncum, klev, klev, icb, nk, inb, t, q, qs, u, &              qs(:ncum, :nl), gz, u(:ncum, :nl), v(:ncum, :nl), p, &
251              v, h, lv, hp, ep, clw, m, sig, ment, qent, uent, vent, nent, &              ph(:ncum, :), th(:ncum, :nl - 1), tv, lv, cpn, ep(:ncum, :), &
252              sij, elij, ments, qents)              clw(:ncum, :), m(:ncum, :), ment(:ncum, :, :), elij(:ncum, :, :), &
253                dtphys, plcl, mp, qp(:ncum, :nl), up(:ncum, :nl), vp(:ncum, :nl), &
254         ! Unsaturated (precipitating) downdrafts              wt(:ncum, :nl), water(:ncum, :nl), evap, b)
255         CALL cv30_unsat(icb(:ncum), inb(:ncum), t, q, qs, gz, u, v, p, ph, th, &         CALL cv30_yield(icb, inb, dtphys, t, q, u, v, gz, p, ph, h, hp, &
256              tv, lv, cpn, ep, sigp, clw, m, ment, elij, delt, plcl, mp, &              lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp(:ncum, 2:nl), &
257              qp(:ncum, :nl), up(:ncum, :nl), vp(:ncum, :nl), wt, water, evap, b)              wt(:ncum, :nl - 1), water(:ncum, :nl), evap, b, ment, qent, uent, &
258                vent, nent, elij, sig, tv, tvp, iflag, precip, VPrecip, ft, fq, &
259         ! Yield (tendencies, precipitation, variables of interface with              fu, fv, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)
        ! other processes, etc)  
        CALL cv30_yield(icb(:ncum), inb(:ncum), delt, t, q, u, v, gz, p, ph, &  
             h, hp, lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp, wt, &  
             water(:ncum, :nl), evap(:ncum, :nl), b, ment, qent, uent, vent, &  
             nent, elij, sig, tv, tvp, iflag, precip, VPrecip, ft, fq, fu, fv, &  
             upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)  
   
260         CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)         CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
261           CALL cv30_uncompress(idcum, iflag, precip, VPrecip, sig, w0, ft, fq, &
262         ! UNCOMPRESS THE FIELDS              fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, da, phi, mp, &
263         iflag1 = 42 ! for non convective points              iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, fu1, fv1, inb1, &
264         CALL cv30_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &              Ma1, upwd1, dnwd1, dnwd01, qcondc1, cape1, da1, phi1, mp1)
             ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, &  
             da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &  
             fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, cape1, da1, &  
             phi1, mp1)  
265      ENDIF      ENDIF
266    
267    end SUBROUTINE cv_driver    end SUBROUTINE cv_driver

Legend:
Removed from v.190  
changed lines
  Added in v.201

  ViewVC Help
Powered by ViewVC 1.1.21