/[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 201 - (show annotations)
Mon Jun 6 17:42:15 2016 UTC (7 years, 11 months ago) by guez
File size: 9996 byte(s)
Removed intermediary objects of cv_thermo_m, access suphec_m
directly. Procedure cv_thermo disappeared, all objects are named
constants.

In cv_driver and below, limited extents of arrays to what is needed.

lv, cpn and th in cv30_compress were set at level nl + 1 but lv1, cpn1
and th1 are not defined at this level. This did not lead to an error
because values at nl + 1 were not used.

Removed test on ok_sync in phystokenc because it is not read at run
time. Printing min and max of output NetCDF variables is heavy and
archaic.

Used histwrite_phy in phytrac.

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 dimphy, ONLY: klev, klon
32
33 real, intent(in):: t1(klon, klev) ! temperature, in K
34 real, intent(in):: q1(klon, klev) ! specific humidity
35 real, intent(in):: qs1(klon, klev) ! saturation specific humidity
36
37 real, intent(in):: u1(klon, klev), v1(klon, klev)
38 ! zonal wind and meridional velocity (m/s)
39
40 real, intent(in):: p1(klon, klev) ! full level pressure, in hPa
41
42 real, intent(in):: ph1(klon, klev + 1)
43 ! Half level pressure, in hPa. These pressures are defined at levels
44 ! 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)
46 ! the first value of the array P1.
47
48 integer, intent(out):: iflag1(:) ! (klon)
49 ! Flag for Emanuel conditions.
50
51 ! 0: Moist convection occurs.
52
53 ! 1: Moist convection occurs, but a CFL condition on the
54 ! subsidence warming is violated. This does not cause the scheme
55 ! to terminate.
56
57 ! 2: Moist convection, but no precipitation because ep(inb) < 1e-4
58
59 ! 3: No moist convection because new cbmf is 0 and old cbmf is 0.
60
61 ! 4: No moist convection; atmosphere is not unstable.
62
63 ! 6: No moist convection because ihmin <= minorig.
64
65 ! 7: No moist convection because unreasonable parcel level
66 ! temperature or specific humidity.
67
68 ! 8: No moist convection: lifted condensation level is above the
69 ! 200 mbar level.
70
71 ! 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)
74 real, intent(out):: fq1(klon, klev) ! specific humidity tendency (s-1)
75
76 real, intent(out):: fu1(klon, klev), fv1(klon, klev)
77 ! forcing (tendency) of zonal and meridional velocity (m/s^2)
78
79 real, intent(out):: precip1(klon) ! convective precipitation rate (mm/day)
80
81 real, intent(out):: VPrecip1(klon, klev + 1)
82 ! vertical profile of convective precipitation (kg/m2/s)
83
84 real, intent(inout):: sig1(klon, klev) ! section of adiabatic updraft
85
86 real, intent(inout):: w01(klon, klev)
87 ! vertical velocity within adiabatic updraft
88
89 integer, intent(out):: icb1(klon)
90 integer, intent(inout):: inb1(klon)
91 real, intent(out):: Ma1(klon, klev) ! mass flux of adiabatic updraft
92
93 real, intent(out):: upwd1(klon, klev)
94 ! total upward mass flux (adiabatic + mixed)
95
96 real, intent(out):: dnwd1(klon, klev) ! saturated downward mass flux (mixed)
97 real, intent(out):: dnwd01(klon, klev) ! unsaturated downward mass flux
98
99 real, intent(out):: qcondc1(klon, klev)
100 ! in-cloud mixing ratio of condensed water
101
102 real, intent(out):: cape1(klon)
103 real, intent(inout):: da1(klon, klev), phi1(klon, klev, klev)
104 real, intent(inout):: mp1(klon, klev)
105
106 ! Local:
107
108 real da(klon, klev), phi(klon, klev, klev)
109 real, allocatable:: mp(:, :) ! (ncum, nl)
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
119 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)
126 real gz1(klon, klev)
127 real hm1(klon, klev)
128 real h1(klon, klev)
129 real tp1(klon, klev)
130 real tvp1(klon, klev)
131 real clw1(klon, klev)
132 real th1(klon, nl) ! potential temperature, in K
133 integer ncum
134
135 ! Compressed fields:
136 integer, allocatable:: idcum(:), iflag(:) ! (ncum)
137 integer, allocatable:: icb(:) ! (ncum)
138 integer nent(klon, klev)
139 integer icbs(klon)
140
141 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)
148 real u(klon, klev), v(klon, klev)
149 real gz(klon, klev), h(klon, klev)
150
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)
160 real pbase(klon), buoybase(klon)
161 real, allocatable:: th(:, :) ! (ncum, nl)
162 real tvp(klon, klev)
163 real sig(klon, klev), w0(klon, klev)
164 real hp(klon, klev), ep(klon, klev)
165 real buoy(klon, klev)
166 real cape(klon)
167 real m(klon, klev), ment(klon, klev, klev), qent(klon, klev, klev)
168 real uent(klon, klev, klev), vent(klon, klev, klev)
169 real ments(klon, klev, klev), qents(klon, klev, klev)
170 real sij(klon, klev, klev), elij(klon, klev, klev)
171 real qp(klon, klev), up(klon, klev), vp(klon, klev)
172 real wt(klon, klev), water(klon, klev)
173 real, allocatable:: evap(:, :) ! (ncum, nl)
174 real, allocatable:: b(:, :) ! (ncum, nl - 1)
175 real ft(klon, klev), fq(klon, klev)
176 real fu(klon, klev), fv(klon, klev)
177 real upwd(klon, klev), dnwd(klon, klev), dnwd0(klon, klev)
178 real Ma(klon, klev), mike(klon, klev), tls(klon, klev)
179 real tps(klon, klev)
180 real precip(klon)
181 real VPrecip(klon, klev + 1)
182 real qcondc(klon, klev) ! cld
183
184 !-------------------------------------------------------------------
185
186 ! SET CONSTANTS AND PARAMETERS
187 CALL cv30_param
188
189 ! INITIALIZE OUTPUT ARRAYS AND PARAMETERS
190
191 do k = 1, klev
192 do i = 1, klon
193 ft1(i, k) = 0.
194 fq1(i, k) = 0.
195 fu1(i, k) = 0.
196 fv1(i, k) = 0.
197 tvp1(i, k) = 0.
198 tp1(i, k) = 0.
199 clw1(i, k) = 0.
200 clw(i, k) = 0.
201 gz1(i, k) = 0.
202 VPrecip1(i, k) = 0.
203 Ma1(i, k) = 0.
204 upwd1(i, k) = 0.
205 dnwd1(i, k) = 0.
206 dnwd01(i, k) = 0.
207 qcondc1(i, k) = 0.
208 end do
209 end do
210
211 precip1 = 0.
212 cape1 = 0.
213 VPrecip1(:, klev + 1) = 0.
214
215 do il = 1, klon
216 sig1(il, klev) = sig1(il, klev) + 1.
217 sig1(il, klev) = min(sig1(il, klev), 12.1)
218 enddo
219
220 CALL cv30_prelim(t1, q1, p1, ph1, lv1, cpn1, tv1, gz1, h1, hm1, th1)
221 CALL cv30_feed(t1, q1, qs1, p1, ph1, gz1, icb1, iflag1, tnk1, qnk1, &
222 gznk1, plcl1)
223 CALL cv30_undilute1(t1, q1, qs1, gz1, plcl1, p1, icb1, tp1, tvp1, clw1, &
224 icbs1)
225 CALL cv30_trigger(icb1, plcl1, p1, th1, tv1, tvp1, pbase1, buoybase1, &
226 iflag1, sig1, w01)
227
228 ncum = count(iflag1 == 0)
229
230 IF (ncum > 0) THEN
231 ! Moist convective adjustment is necessary
232 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, &
238 cpn1, p1, ph1, tv1, tp1, tvp1, clw1, sig1, w01, icb, icbs, plcl, &
239 tnk, qnk, gznk, pbase, buoybase, t, q, qs, u, v, gz, th, h, lv, &
240 cpn, p, ph, tv, tp, tvp, clw, sig, w0)
241 CALL cv30_undilute2(icb, icbs(:ncum), tnk, qnk, gznk, t, qs, gz, p, h, &
242 tv, lv, pbase(:ncum), buoybase(:ncum), plcl, inb, tp, tvp, &
243 clw, hp, ep, buoy)
244 CALL cv30_closure(icb, inb, pbase, p, ph(:ncum, :), tv, buoy, &
245 sig, w0, cape, m)
246 CALL cv30_mixing(icb, inb, t, q, qs, u, v, h, lv, &
247 hp, ep, clw, m, sig, ment, qent, uent, vent, nent, sij, elij, &
248 ments, qents)
249 CALL cv30_unsat(icb, inb, t(:ncum, :nl), q(:ncum, :nl), &
250 qs(:ncum, :nl), gz, u(:ncum, :nl), v(:ncum, :nl), p, &
251 ph(:ncum, :), th(:ncum, :nl - 1), tv, lv, cpn, ep(:ncum, :), &
252 clw(:ncum, :), m(:ncum, :), ment(:ncum, :, :), elij(:ncum, :, :), &
253 dtphys, plcl, mp, qp(:ncum, :nl), up(:ncum, :nl), vp(:ncum, :nl), &
254 wt(:ncum, :nl), water(:ncum, :nl), evap, b)
255 CALL cv30_yield(icb, inb, dtphys, t, q, u, v, gz, p, ph, h, hp, &
256 lv, cpn, th, ep, clw, m, tp, mp, qp, up, vp(:ncum, 2:nl), &
257 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 fu, fv, upwd, dnwd, dnwd0, ma, mike, tls, tps, qcondc)
260 CALL cv30_tracer(klon, ncum, klev, ment, sij, da, phi)
261 CALL cv30_uncompress(idcum, iflag, precip, VPrecip, sig, w0, ft, fq, &
262 fu, fv, inb, Ma, upwd, dnwd, dnwd0, qcondc, cape, da, phi, mp, &
263 iflag1, precip1, VPrecip1, sig1, w01, ft1, fq1, fu1, fv1, inb1, &
264 Ma1, upwd1, dnwd1, dnwd01, qcondc1, cape1, da1, phi1, mp1)
265 ENDIF
266
267 end SUBROUTINE cv_driver
268
269 end module cv_driver_m

  ViewVC Help
Powered by ViewVC 1.1.21