/[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 189 - (hide annotations)
Tue Mar 29 15:20:23 2016 UTC (8 years, 2 months ago) by guez
File size: 9954 byte(s)
There was a function gr_phy_write_3d in dyn3d and a function
gr_phy_write_2d in module grid_change. Moved them into a new module
gr_phy_write_m under a generic interface gr_phy_write. Replaced calls
to gr_fi_ecrit by calls to gr_phy_write.

Removed arguments len, nloc and nd of cv30_compress.

Removed arguments wd and wd1 of cv30_uncompress, wd of cv30_yield, wd
of concvl, wd1 of cv_driver. Was just filled with 0. Removed option
ok_gust in physiq, never used.

In cv30_unsat, cv30_yield and cv_driver, we only need to define b to
level nl - 1.

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

  ViewVC Help
Powered by ViewVC 1.1.21