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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 213 - (show annotations)
Mon Feb 27 15:44:55 2017 UTC (7 years, 2 months ago) by guez
File size: 5905 byte(s)
Removed module conema3_m. Moved variables epmax and iflag_clw of
conema3_m to conf_phys_m, where they are defined. Removed unused
variable ok_adj_ema of conema3_m.

Added variables d_t_ec, dtsw0 and dtlw0 to histins.nc (following LMDZ).

Removed case not lessivage in phytrac. (Not used in LMDZ without INCA
either.)

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

  ViewVC Help
Powered by ViewVC 1.1.21