/[lmdze]/trunk/Sources/phylmd/CV30_routines/cv30_undilute2.f
ViewVC logotype

Annotation of /trunk/Sources/phylmd/CV30_routines/cv30_undilute2.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 196 - (hide annotations)
Mon May 23 13:50:39 2016 UTC (8 years ago) by guez
File size: 5864 byte(s)
Removed argument icbmax of cv30_feed, not used in cv_driver (not used
in LMDZ either).

Clearer to have iflag1 = 0 in cv30_feed than in cv_driver. Clearer to
have iflag1 = 42 in cv30_uncompress than in cv_driver.

Removed argument iflag of cv30_compress. Since there is iflag1 = 42 in
cv30_uncompress, iflag1 that comes out of cv_driver is disconnected
from iflag1 computed by cv30_feed and cv30_trigger.

Removed some references to convect3 and convect4 in comments. This
program derives from the convect3 version, we do not need to know
about other versions.

Bug fix in cv30_undilute1: icbs1 was not made >= 2.

1 guez 185 module cv30_undilute2_m
2 guez 47
3 guez 183 implicit none
4 guez 47
5 guez 183 contains
6 guez 47
7 guez 195 SUBROUTINE cv30_undilute2(icb, icbs, nk, tnk, qnk, gznk, t, qs, gz, p, h, &
8     tv, lv, pbase, buoybase, plcl, inb, tp, tvp, clw, hp, ep, buoy)
9 guez 47
10 guez 187 ! Undilute (adiabatic) updraft, second part. Purpose: find the
11     ! rest of the lifted parcel temperatures; compute the
12     ! precipitation efficiencies and the fraction of precipitation
13     ! falling outside of cloud; find the level of neutral buoyancy.
14 guez 47
15 guez 183 ! Vertical profile of buoyancy computed here (use of buoybase).
16 guez 47
17 guez 183 use conema3_m, only: epmax
18 guez 195 use cv30_param_m, only: minorig, nl
19 guez 190 use cv_thermo_m, only: cl, clmcpv, cpd, cpv, eps, lv0, rrv
20 guez 187 USE dimphy, ONLY: klon, klev
21 guez 47
22 guez 195 integer, intent(in):: icb(:), icbs(:) ! (ncum)
23 guez 187 ! icbs is the first level above LCL (may differ from icb)
24    
25     integer, intent(in):: nk(klon)
26     real, intent(in):: tnk(klon), qnk(klon), gznk(klon)
27     real, intent(in):: t(klon, klev), qs(klon, klev), gz(klon, klev)
28     real, intent(in):: p(klon, klev), h(klon, klev)
29     real, intent(in):: tv(klon, klev), lv(klon, klev)
30 guez 196 real, intent(in):: pbase(:), buoybase(:), plcl(:) ! (ncum)
31 guez 187
32 guez 183 ! outputs:
33 guez 187 integer, intent(out):: inb(:) ! (ncum)
34     ! first model level above the level of neutral buoyancy of the
35 guez 192 ! parcel (1 <= inb <= nl - 1)
36 guez 187
37     real tp(klon, klev), tvp(klon, klev), clw(klon, klev)
38 guez 183 ! condensed water not removed from tvp
39 guez 195 real hp(klon, klev), ep(klon, klev)
40 guez 187 real buoy(klon, klev)
41 guez 47
42 guez 183 ! Local:
43 guez 195
44     integer ncum
45    
46     real, parameter:: pbcrit = 150.
47     ! critical cloud depth (mbar) beneath which the precipitation
48     ! efficiency is assumed to be zero
49    
50     real, parameter:: ptcrit = 500.
51     ! cloud depth (mbar) above which the precipitation efficiency is
52     ! assumed to be unity
53    
54     real, parameter:: dtovsh = - 0.2 ! dT for overshoot
55    
56 guez 183 integer i, k
57     real tg, qg, ahg, alv, s, tc, es, denom
58     real pden
59 guez 187 real ah0(klon)
60 guez 47
61 guez 183 !---------------------------------------------------------------------
62 guez 47
63 guez 195 ncum = size(icb)
64    
65 guez 183 ! SOME INITIALIZATIONS
66 guez 47
67 guez 183 do k = 1, nl
68     do i = 1, ncum
69 guez 195 ep(i, k) = 0.
70 guez 183 end do
71     end do
72 guez 47
73 guez 183 ! FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
74 guez 47
75 guez 183 ! The procedure is to solve the equation.
76     ! cp * tp + L * qp + phi = cp * tnk + L * qnk + gznk.
77 guez 47
78 guez 183 ! Calculate certain parcel quantities, including static energy
79 guez 47
80 guez 183 do i = 1, ncum
81     ah0(i) = (cpd * (1. - qnk(i)) + cl * qnk(i)) * tnk(i) &
82     + qnk(i) * (lv0 - clmcpv * (tnk(i) - 273.15)) + gznk(i)
83     end do
84 guez 47
85 guez 183 ! Find lifted parcel quantities above cloud base
86 guez 47
87 guez 183 do k = minorig + 1, nl
88     do i = 1, ncum
89     if (k >= (icbs(i) + 1)) then
90     tg = t(i, k)
91     qg = qs(i, k)
92     alv = lv0 - clmcpv * (t(i, k) - 273.15)
93 guez 47
94 guez 183 ! First iteration.
95 guez 47
96 guez 183 s = cpd * (1. - qnk(i)) + cl * qnk(i) &
97     + alv * alv * qg / (rrv * t(i, k) * t(i, k))
98     s = 1. / s
99 guez 47
100 guez 183 ahg = cpd * tg + (cl - cpd) * qnk(i) * tg + alv * qg + gz(i, k)
101     tg = tg + s * (ah0(i) - ahg)
102 guez 47
103 guez 183 tc = tg - 273.15
104     denom = 243.5 + tc
105 guez 195 denom = MAX(denom, 1.)
106 guez 47
107 guez 183 es = 6.112 * exp(17.67 * tc / denom)
108 guez 47
109 guez 183 qg = eps * es / (p(i, k) - es * (1. - eps))
110 guez 47
111 guez 183 ! Second iteration.
112 guez 47
113 guez 183 ahg = cpd * tg + (cl - cpd) * qnk(i) * tg + alv * qg + gz(i, k)
114     tg = tg + s * (ah0(i) - ahg)
115 guez 47
116 guez 183 tc = tg - 273.15
117     denom = 243.5 + tc
118 guez 195 denom = MAX(denom, 1.)
119 guez 47
120 guez 183 es = 6.112 * exp(17.67 * tc / denom)
121    
122     qg = eps * es / (p(i, k) - es * (1. - eps))
123    
124     alv = lv0 - clmcpv * (t(i, k) - 273.15)
125    
126     ! no approximation:
127     tp(i, k) = (ah0(i) - gz(i, k) - alv * qg) &
128     / (cpd + (cl - cpd) * qnk(i))
129    
130     clw(i, k) = qnk(i) - qg
131 guez 195 clw(i, k) = max(0., clw(i, k))
132 guez 183 ! qg utilise au lieu du vrai mixing ratio rg:
133     tvp(i, k) = tp(i, k) * (1. + qg / eps - qnk(i)) ! whole thing
134     endif
135     end do
136     end do
137    
138 guez 195 ! SET THE PRECIPITATION EFFICIENCIES
139     ! It MAY BE a FUNCTION OF TP(I), P(I) AND CLW(I)
140 guez 183 do k = 1, nl
141     do i = 1, ncum
142     pden = ptcrit - pbcrit
143     ep(i, k) = (plcl(i) - p(i, k) - pbcrit) / pden * epmax
144 guez 195 ep(i, k) = max(ep(i, k), 0.)
145 guez 192 ep(i, k) = min(ep(i, k), epmax)
146 guez 183 end do
147     end do
148    
149     ! CALCULATE VIRTUAL TEMPERATURE AND LIFTED PARCEL
150     ! VIRTUAL TEMPERATURE
151    
152     ! tvp est calcule en une seule fois, et sans retirer
153     ! l'eau condensee (~> reversible CAPE)
154     do i = 1, ncum
155 guez 186 tp(i, nl + 1) = tp(i, nl)
156 guez 183 end do
157    
158     ! EFFECTIVE VERTICAL PROFILE OF BUOYANCY:
159    
160     ! first estimate of buoyancy:
161     do i = 1, ncum
162     do k = 1, nl
163     buoy(i, k) = tvp(i, k) - tv(i, k)
164     end do
165     end do
166    
167     ! set buoyancy = buoybase for all levels below base
168     ! for safety, set buoy(icb) = buoybase
169     do i = 1, ncum
170     do k = 1, nl
171     if ((k >= icb(i)) .and. (k <= nl) .and. (p(i, k) >= pbase(i))) then
172     buoy(i, k) = buoybase(i)
173     endif
174     end do
175     buoy(icb(i), k) = buoybase(i)
176     end do
177    
178 guez 187 ! Compute inb:
179 guez 183
180 guez 187 inb = nl - 1
181 guez 183
182     do i = 1, ncum
183     do k = 1, nl - 1
184     if ((k >= icb(i)) .and. (buoy(i, k) < dtovsh)) then
185     inb(i) = MIN(inb(i), k)
186     endif
187     end do
188     end do
189    
190     ! CALCULATE LIQUID WATER STATIC ENERGY OF LIFTED PARCEL
191    
192 guez 186 do k = 1, nl + 1
193 guez 183 do i = 1, ncum
194     hp(i, k) = h(i, k)
195     enddo
196     enddo
197    
198     do k = minorig + 1, nl
199     do i = 1, ncum
200     if (k >= icb(i) .and. k <= inb(i)) hp(i, k) = h(i, nk(i)) &
201     + (lv(i, k) + (cpd - cpv) * t(i, k)) * ep(i, k) * clw(i, k)
202     end do
203     end do
204    
205 guez 185 end SUBROUTINE cv30_undilute2
206 guez 183
207 guez 185 end module cv30_undilute2_m

  ViewVC Help
Powered by ViewVC 1.1.21