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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 195 - (hide annotations)
Wed May 18 17:56:44 2016 UTC (8 years ago) by guez
File size: 9626 byte(s)
In cv30_feed, iflag1 is 0 on entry so we can simplify the test for
iflag1 = 7.

In cv30_feed, for the computation of icb, replaced sequential search
(with a useless end of loop on k) by a call to locate.

In CV30 routines, replaced len, nloc, nd, na by klon or
klev. Philosophy: no more generality than actually necessary.

Converted as many variables as possible to named constants in
cv30_param_m and downgraded pbcrit, ptcrit, dtovsh, dpbase, dttrig,
tau, delta to local objects in procedures. spfac, betad and omtrain
are useless and removed.

Instead of filling the array sigp with the constant spfac in
cv30_undilute2, just made sigp a constant in cv30_unsat.

In cv_driver, define as allocatable variables that are only
used on the range (ncum, nl).

1 guez 52 module cv_driver_m
2 guez 3
3 guez 52 implicit none
4 guez 3
5 guez 52 contains
6 guez 3
7 guez 183 SUBROUTINE cv_driver(t1, q1, qs1, u1, v1, p1, ph1, iflag1, ft1, fq1, fu1, &
8 guez 195 fv1, precip1, VPrecip1, sig1, w01, icb1, inb1, Ma1, upwd1, dnwd1, &
9 guez 189 dnwd01, qcondc1, cape1, da1, phi1, mp1)
10 guez 3
11 guez 91 ! From LMDZ4/libf/phylmd/cv_driver.F, version 1.3, 2005/04/15 12:36:17
12 guez 69 ! Main driver for convection
13 guez 97 ! Author: S. Bony, March 2002
14 guez 69
15 guez 72 ! Several modules corresponding to different physical processes
16    
17 guez 195 use comconst, only: dtphys
18 guez 187 use cv30_closure_m, only: cv30_closure
19 guez 185 use cv30_compress_m, only: cv30_compress
20     use cv30_feed_m, only: cv30_feed
21     use cv30_mixing_m, only: cv30_mixing
22 guez 188 use cv30_param_m, only: cv30_param, nl
23 guez 185 use cv30_prelim_m, only: cv30_prelim
24     use cv30_tracer_m, only: cv30_tracer
25 guez 189 use cv30_trigger_m, only: cv30_trigger
26 guez 185 use cv30_uncompress_m, only: cv30_uncompress
27 guez 195 use cv30_undilute1_m, only: cv30_undilute1
28 guez 185 use cv30_undilute2_m, only: cv30_undilute2
29     use cv30_unsat_m, only: cv30_unsat
30     use cv30_yield_m, only: cv30_yield
31 guez 190 use cv_thermo_m, only: cv_thermo
32 guez 62 USE dimphy, ONLY: klev, klon
33    
34 guez 189 real, intent(in):: t1(klon, klev) ! temperature (K)
35     real, intent(in):: q1(klon, klev) ! specific humidity
36     real, intent(in):: qs1(klon, klev) ! saturation specific humidity
37 guez 103
38 guez 187 real, intent(in):: u1(klon, klev), v1(klon, klev)
39 guez 189 ! zonal wind and meridional velocity (m/s)
40 guez 72
41 guez 189 real, intent(in):: p1(klon, klev) ! full level pressure (hPa)
42 guez 180
43 guez 187 real, intent(in):: ph1(klon, klev + 1)
44 guez 189 ! Half level pressure (hPa). These pressures are defined at levels
45     ! intermediate between those of P1, T1, Q1 and QS1. The first
46     ! value of PH should be greater than (i.e. at a lower level than)
47     ! the first value of the array P1.
48 guez 180
49 guez 187 integer, intent(out):: iflag1(klon)
50     ! Flag for Emanuel conditions.
51 guez 3
52 guez 187 ! 0: Moist convection occurs.
53 guez 3
54 guez 187 ! 1: Moist convection occurs, but a CFL condition on the
55     ! subsidence warming is violated. This does not cause the scheme
56     ! to terminate.
57 guez 3
58 guez 187 ! 2: Moist convection, but no precipitation because ep(inb) < 1e-4
59 guez 103
60 guez 187 ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
61 guez 62
62 guez 195 ! 4: No moist convection; atmosphere is not unstable.
63 guez 62
64 guez 195 ! 6: No moist convection because ihmin <= minorig.
65 guez 62
66 guez 187 ! 7: No moist convection because unreasonable parcel level
67     ! temperature or specific humidity.
68 guez 62
69 guez 187 ! 8: No moist convection: lifted condensation level is above the
70 guez 195 ! 200 mbar level.
71 guez 62
72 guez 195 ! 9: No moist convection: cloud base is higher than the level NL-1.
73 guez 62
74 guez 189 real, intent(out):: ft1(klon, klev) ! temperature tendency (K/s)
75     real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1)
76 guez 62
77 guez 187 real, intent(out):: fu1(klon, klev), fv1(klon, klev)
78 guez 189 ! forcing (tendency) of zonal and meridional velocity (m/s^2)
79 guez 62
80 guez 187 real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day)
81 guez 62
82 guez 187 real, intent(out):: VPrecip1(klon, klev + 1)
83     ! vertical profile of convective precipitation (kg/m2/s)
84 guez 62
85 guez 189 real, intent(inout):: sig1(klon, klev) ! section of adiabatic updraft
86 guez 62
87 guez 187 real, intent(inout):: w01(klon, klev)
88     ! vertical velocity within adiabatic updraft
89 guez 62
90 guez 187 integer, intent(out):: icb1(klon)
91     integer, intent(inout):: inb1(klon)
92 guez 189 real, intent(out):: Ma1(klon, klev) ! mass flux of adiabatic updraft
93 guez 62
94 guez 187 real, intent(out):: upwd1(klon, klev)
95 guez 189 ! total upward mass flux (adiabatic + mixed)
96 guez 62
97 guez 187 real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
98     real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
99 guez 62
100 guez 189 real, intent(out):: qcondc1(klon, klev)
101     ! in-cloud mixing ratio of condensed water
102 guez 62
103 guez 189 real, intent(out):: cape1(klon)
104 guez 187 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
105     real, intent(inout):: mp1(klon, klev)
106 guez 62
107 guez 187 ! Local:
108 guez 62
109 guez 103 real da(klon, klev), phi(klon, klev, klev), mp(klon, klev)
110 guez 97 integer i, k, il
111 guez 52 integer icbmax
112     integer nk1(klon)
113     integer icbs1(klon)
114     real plcl1(klon)
115     real tnk1(klon)
116     real qnk1(klon)
117     real gznk1(klon)
118     real pbase1(klon)
119     real buoybase1(klon)
120     real lv1(klon, klev)
121     real cpn1(klon, klev)
122     real tv1(klon, klev)
123     real gz1(klon, klev)
124     real hm1(klon, klev)
125     real h1(klon, klev)
126     real tp1(klon, klev)
127     real tvp1(klon, klev)
128     real clw1(klon, klev)
129     real th1(klon, klev)
130     integer ncum
131 guez 62
132 guez 187 ! Compressed fields:
133 guez 103 integer idcum(klon)
134 guez 195 integer iflag(klon), nk(klon)
135     integer, allocatable:: icb(:) ! (ncum)
136 guez 103 integer nent(klon, klev)
137     integer icbs(klon)
138 guez 183 integer inb(klon)
139 guez 182 real plcl(klon), tnk(klon), qnk(klon), gznk(klon)
140 guez 103 real t(klon, klev), q(klon, klev), qs(klon, klev)
141     real u(klon, klev), v(klon, klev)
142     real gz(klon, klev), h(klon, klev), lv(klon, klev), cpn(klon, klev)
143 guez 195 real p(klon, klev) ! pressure at full level, in hPa
144     real ph(klon, klev + 1), tv(klon, klev), tp(klon, klev)
145 guez 103 real clw(klon, klev)
146     real pbase(klon), buoybase(klon), th(klon, klev)
147     real tvp(klon, klev)
148     real sig(klon, klev), w0(klon, klev)
149 guez 195 real hp(klon, klev), ep(klon, klev)
150 guez 183 real buoy(klon, klev)
151 guez 103 real cape(klon)
152     real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
153     real uent(klon, klev, klev), vent(klon, klev, klev)
154     real ments(klon, klev, klev), qents(klon, klev, klev)
155     real sij(klon, klev, klev), elij(klon, klev, klev)
156     real qp(klon, klev), up(klon, klev), vp(klon, klev)
157 guez 195 real wt(klon, klev), water(klon, klev)
158     real, allocatable:: evap(:, :) ! (ncum, nl)
159 guez 189 real, allocatable:: b(:, :) ! (ncum, nl - 1)
160 guez 188 real ft(klon, klev), fq(klon, klev)
161 guez 103 real fu(klon, klev), fv(klon, klev)
162     real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
163     real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
164 guez 183 real tps(klon, klev)
165 guez 103 real precip(klon)
166 guez 180 real VPrecip(klon, klev + 1)
167 guez 103 real qcondc(klon, klev) ! cld
168 guez 3
169 guez 52 !-------------------------------------------------------------------
170 guez 3
171 guez 180 ! SET CONSTANTS AND PARAMETERS
172    
173     ! set thermodynamical constants:
174 guez 69 CALL cv_thermo
175 guez 52
176 guez 195 CALL cv30_param
177 guez 52
178 guez 180 ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
179 guez 3
180 guez 103 do k = 1, klev
181     do i = 1, klon
182 guez 187 ft1(i, k) = 0.
183     fq1(i, k) = 0.
184     fu1(i, k) = 0.
185     fv1(i, k) = 0.
186     tvp1(i, k) = 0.
187     tp1(i, k) = 0.
188     clw1(i, k) = 0.
189     clw(i, k) = 0.
190 guez 103 gz1(i, k) = 0.
191 guez 52 VPrecip1(i, k) = 0.
192 guez 187 Ma1(i, k) = 0.
193     upwd1(i, k) = 0.
194     dnwd1(i, k) = 0.
195     dnwd01(i, k) = 0.
196     qcondc1(i, k) = 0.
197 guez 52 end do
198     end do
199 guez 3
200 guez 195 precip1 = 0.
201     iflag1 = 0
202     cape1 = 0.
203     VPrecip1(:, klev + 1) = 0.
204 guez 3
205 guez 181 do il = 1, klon
206     sig1(il, klev) = sig1(il, klev) + 1.
207     sig1(il, klev) = min(sig1(il, klev), 12.1)
208     enddo
209 guez 3
210 guez 185 CALL cv30_prelim(klon, klev, klev + 1, t1, q1, p1, ph1, lv1, cpn1, tv1, &
211 guez 181 gz1, h1, hm1, th1)
212 guez 195 CALL cv30_feed(t1, q1, qs1, p1, ph1, gz1, nk1, icb1, icbmax, iflag1, &
213     tnk1, qnk1, gznk1, plcl1)
214     CALL cv30_undilute1(t1, q1, qs1, gz1, plcl1, p1, nk1, icb1, tp1, tvp1, &
215     clw1, icbs1)
216     CALL cv30_trigger(icb1, plcl1, p1, th1, tv1, tvp1, pbase1, buoybase1, &
217     iflag1, sig1, w01)
218 guez 3
219 guez 180 ! Moist convective adjustment is necessary
220 guez 3
221 guez 91 ncum = 0
222 guez 103 do i = 1, klon
223 guez 180 if (iflag1(i) == 0) then
224     ncum = ncum + 1
225 guez 91 idcum(ncum) = i
226 guez 52 endif
227     end do
228 guez 3
229 guez 103 IF (ncum > 0) THEN
230 guez 195 allocate(b(ncum, nl - 1), evap(ncum, nl), icb(ncum))
231 guez 189 CALL cv30_compress(ncum, iflag1, nk1, icb1, icbs1, plcl1, tnk1, qnk1, &
232     gznk1, pbase1, buoybase1, t1, q1, qs1, u1, v1, gz1, th1, h1, lv1, &
233     cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, iflag, nk, icb, &
234     icbs, plcl, tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, &
235     th, h, lv, cpn, p, ph, tv, tp, tvp, clw, sig, w0)
236 guez 195 CALL cv30_undilute2(icb, icbs(:ncum), nk, tnk, qnk, gznk, t, qs, gz, &
237     p, h, tv, lv, pbase, buoybase, plcl, inb(:ncum), tp, tvp, clw, &
238     hp, ep, buoy)
239     CALL cv30_closure(icb, inb(:ncum), pbase, p, ph, tv, buoy, sig, w0, &
240     cape, m)
241     CALL cv30_mixing(icb, nk(:ncum), inb(:ncum), t, q, qs, u, v, h, lv, &
242     hp, ep, clw, m, sig, ment, qent, uent, vent, nent, sij, elij, &
243     ments, qents)
244     CALL cv30_unsat(icb, inb(:ncum), t(:ncum, :nl), q(:ncum, :nl), &
245     qs(:ncum, :nl), gz, u, v, p, ph, th(:ncum, :nl - 1), tv, lv, cpn, &
246     ep(:ncum, :), clw(:ncum, :), m(:ncum, :), ment(:ncum, :, :), &
247     elij(:ncum, :, :), dtphys, plcl, mp, qp(:ncum, :nl), &
248     up(:ncum, :nl), vp(:ncum, :nl), wt(:ncum, :nl), &
249     water(:ncum, :nl), evap, b)
250     CALL cv30_yield(icb, inb(:ncum), dtphys, t, q, u, v, gz, p, ph, h, hp, &
251     lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp(:ncum, 2:nl), &
252     wt(:ncum, :nl - 1), water(:ncum, :nl), evap, b, ment, qent, uent, &
253     vent, nent, elij, sig, tv, tvp, iflag, precip, VPrecip, ft, fq, &
254     fu, fv, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)
255 guez 185 CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
256 guez 3
257 guez 180 ! UNCOMPRESS THE FIELDS
258 guez 189 iflag1 = 42 ! for non convective points
259 guez 185 CALL cv30_uncompress(idcum(:ncum), iflag, precip, VPrecip, sig, w0, &
260 guez 189 ft, fq, fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, &
261 guez 181 da, phi, mp, iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, &
262 guez 189 fu1, fv1, inb1, Ma1, upwd1, dnwd1, dnwd01, qcondc1, cape1, da1, &
263     phi1, mp1)
264 guez 183 ENDIF
265 guez 3
266 guez 52 end SUBROUTINE cv_driver
267 guez 3
268 guez 52 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21