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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 198 - (show annotations)
Tue May 31 16:17:35 2016 UTC (7 years, 11 months ago) by guez
File size: 9456 byte(s)
Removed variables nk1 and nk in cv_driver and below. These arrays were
just equal to the constant minorig. (This is also the case in LMDZ.)

In cv_thermo, removed some variables which were copies of variables of
suphec_m. Changed some variables to named constants.

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

  ViewVC Help
Powered by ViewVC 1.1.21