/[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 150 - (hide annotations)
Thu Jun 18 13:49:26 2015 UTC (8 years, 11 months ago) by guez
Original Path: trunk/Sources/phylmd/CV3_routines/cv3_undilute2.f
File size: 11450 byte(s)
Removed unused arguments of groupe, cv3_undilute2, cv_undilute2,
interfsur_lim, drag_noro, orodrag, gwprofil

Chickened out of revision 148: back to double precision in
invert_zoom_x (and overloaded rtsafe).

1 guez 47
2     SUBROUTINE cv3_undilute2(nloc,ncum,nd,icb,icbs,nk &
3 guez 150 ,tnk,qnk,gznk,t,qs,gz &
4 guez 47 ,p,h,tv,lv,pbase,buoybase,plcl &
5     ,inb,tp,tvp,clw,hp,ep,sigp,buoy)
6     use conema3_m
7 guez 69 use cv3_param_m
8 guez 47 use cvthermo
9     implicit none
10    
11     !---------------------------------------------------------------------
12     ! Purpose:
13     ! FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
14     ! &
15     ! COMPUTE THE PRECIPITATION EFFICIENCIES AND THE
16     ! FRACTION OF PRECIPITATION FALLING OUTSIDE OF CLOUD
17     ! &
18     ! FIND THE LEVEL OF NEUTRAL BUOYANCY
19     !
20     ! Main differences convect3/convect4:
21     ! - icbs (input) is the first level above LCL (may differ from icb)
22     ! - many minor differences in the iterations
23     ! - condensed water not removed from tvp in convect3
24     ! - vertical profile of buoyancy computed here (use of buoybase)
25     ! - the determination of inb is different
26     ! - no inb1, only inb in output
27     !---------------------------------------------------------------------
28    
29    
30     ! inputs:
31 guez 97 integer, intent(in):: ncum, nd, nloc
32 guez 47 integer icb(nloc), icbs(nloc), nk(nloc)
33 guez 150 real t(nloc,nd), qs(nloc,nd), gz(nloc,nd)
34 guez 47 real p(nloc,nd)
35     real tnk(nloc), qnk(nloc), gznk(nloc)
36     real lv(nloc,nd), tv(nloc,nd), h(nloc,nd)
37     real pbase(nloc), buoybase(nloc), plcl(nloc)
38    
39     ! outputs:
40     integer inb(nloc)
41     real tp(nloc,nd), tvp(nloc,nd), clw(nloc,nd)
42     real ep(nloc,nd), sigp(nloc,nd), hp(nloc,nd)
43     real buoy(nloc,nd)
44    
45     ! local variables:
46     integer i, k
47 guez 150 real tg,qg,ahg,alv,s,tc,es,denom,rg
48     real pden
49     real ah0(nloc)
50 guez 47
51     !=====================================================================
52     ! --- SOME INITIALIZATIONS
53     !=====================================================================
54    
55     do 170 k=1,nl
56     do 160 i=1,ncum
57     ep(i,k)=0.0
58     sigp(i,k)=spfac
59     160 continue
60     170 continue
61    
62     !=====================================================================
63     ! --- FIND THE REST OF THE LIFTED PARCEL TEMPERATURES
64     !=====================================================================
65     !
66     ! --- The procedure is to solve the equation.
67     ! cp*tp+L*qp+phi=cp*tnk+L*qnk+gznk.
68     !
69     ! *** Calculate certain parcel quantities, including static energy ***
70     !
71     !
72     do 240 i=1,ncum
73     ah0(i)=(cpd*(1.-qnk(i))+cl*qnk(i))*tnk(i) &
74     +qnk(i)*(lv0-clmcpv*(tnk(i)-273.15))+gznk(i)
75     240 continue
76     !
77     !
78     ! *** Find lifted parcel quantities above cloud base ***
79     !
80     !
81     do 300 k=minorig+1,nl
82     do 290 i=1,ncum
83     ! ori if(k.ge.(icb(i)+1))then
84     if(k.ge.(icbs(i)+1))then ! convect3
85     tg=t(i,k)
86     qg=qs(i,k)
87     !debug alv=lv0-clmcpv*(t(i,k)-t0)
88     alv=lv0-clmcpv*(t(i,k)-273.15)
89     !
90     ! First iteration.
91     !
92     ! ori s=cpd+alv*alv*qg/(rrv*t(i,k)*t(i,k))
93     s=cpd*(1.-qnk(i))+cl*qnk(i) &
94     +alv*alv*qg/(rrv*t(i,k)*t(i,k)) ! convect3
95     s=1./s
96     ! ori ahg=cpd*tg+(cl-cpd)*qnk(i)*t(i,k)+alv*qg+gz(i,k)
97     ahg=cpd*tg+(cl-cpd)*qnk(i)*tg+alv*qg+gz(i,k) ! convect3
98     tg=tg+s*(ah0(i)-ahg)
99     ! ori tg=max(tg,35.0)
100     !debug tc=tg-t0
101     tc=tg-273.15
102     denom=243.5+tc
103     denom=MAX(denom,1.0) ! convect3
104     ! ori if(tc.ge.0.0)then
105     es=6.112*exp(17.67*tc/denom)
106     ! ori else
107     ! ori es=exp(23.33086-6111.72784/tg+0.15215*log(tg))
108     ! ori endif
109     qg=eps*es/(p(i,k)-es*(1.-eps))
110     !
111     ! Second iteration.
112     !
113     ! ori s=cpd+alv*alv*qg/(rrv*t(i,k)*t(i,k))
114     ! ori s=1./s
115     ! ori ahg=cpd*tg+(cl-cpd)*qnk(i)*t(i,k)+alv*qg+gz(i,k)
116     ahg=cpd*tg+(cl-cpd)*qnk(i)*tg+alv*qg+gz(i,k) ! convect3
117     tg=tg+s*(ah0(i)-ahg)
118     ! ori tg=max(tg,35.0)
119     !debug tc=tg-t0
120     tc=tg-273.15
121     denom=243.5+tc
122     denom=MAX(denom,1.0) ! convect3
123     ! ori if(tc.ge.0.0)then
124     es=6.112*exp(17.67*tc/denom)
125     ! ori else
126     ! ori es=exp(23.33086-6111.72784/tg+0.15215*log(tg))
127     ! ori endif
128     qg=eps*es/(p(i,k)-es*(1.-eps))
129     !
130     !debug alv=lv0-clmcpv*(t(i,k)-t0)
131     alv=lv0-clmcpv*(t(i,k)-273.15)
132     ! print*,'cpd dans convect2 ',cpd
133     ! print*,'tp(i,k),ah0(i),cl,cpd,qnk(i),t(i,k),gz(i,k),alv,qg,cpd'
134     ! print*,tp(i,k),ah0(i),cl,cpd,qnk(i),t(i,k),gz(i,k),alv,qg,cpd
135    
136     ! ori c approximation here:
137     ! ori tp(i,k)=(ah0(i)-(cl-cpd)*qnk(i)*t(i,k)-gz(i,k)-alv*qg)/cpd
138    
139     ! convect3: no approximation:
140     tp(i,k)=(ah0(i)-gz(i,k)-alv*qg)/(cpd+(cl-cpd)*qnk(i))
141    
142     clw(i,k)=qnk(i)-qg
143     clw(i,k)=max(0.0,clw(i,k))
144     rg=qg/(1.-qnk(i))
145     ! ori tvp(i,k)=tp(i,k)*(1.+rg*epsi)
146     ! convect3: (qg utilise au lieu du vrai mixing ratio rg):
147     tvp(i,k)=tp(i,k)*(1.+qg/eps-qnk(i)) ! whole thing
148     endif
149     290 continue
150     300 continue
151     !
152     !=====================================================================
153     ! --- SET THE PRECIPITATION EFFICIENCIES AND THE FRACTION OF
154     ! --- PRECIPITATION FALLING OUTSIDE OF CLOUD
155     ! --- THESE MAY BE FUNCTIONS OF TP(I), P(I) AND CLW(I)
156     !=====================================================================
157     !
158     ! ori do 320 k=minorig+1,nl
159     do 320 k=1,nl ! convect3
160     do 310 i=1,ncum
161     pden=ptcrit-pbcrit
162     ep(i,k)=(plcl(i)-p(i,k)-pbcrit)/pden*epmax
163     ep(i,k)=amax1(ep(i,k),0.0)
164     ep(i,k)=amin1(ep(i,k),epmax)
165     sigp(i,k)=spfac
166     ! ori if(k.ge.(nk(i)+1))then
167     ! ori tca=tp(i,k)-t0
168     ! ori if(tca.ge.0.0)then
169     ! ori elacrit=elcrit
170     ! ori else
171     ! ori elacrit=elcrit*(1.0-tca/tlcrit)
172     ! ori endif
173     ! ori elacrit=max(elacrit,0.0)
174     ! ori ep(i,k)=1.0-elacrit/max(clw(i,k),1.0e-8)
175     ! ori ep(i,k)=max(ep(i,k),0.0 )
176     ! ori ep(i,k)=min(ep(i,k),1.0 )
177     ! ori sigp(i,k)=sigs
178     ! ori endif
179     310 continue
180     320 continue
181     !
182     !=====================================================================
183     ! --- CALCULATE VIRTUAL TEMPERATURE AND LIFTED PARCEL
184     ! --- VIRTUAL TEMPERATURE
185     !=====================================================================
186     !
187     ! dans convect3, tvp est calcule en une seule fois, et sans retirer
188     ! l'eau condensee (~> reversible CAPE)
189     !
190     ! ori do 340 k=minorig+1,nl
191     ! ori do 330 i=1,ncum
192     ! ori if(k.ge.(icb(i)+1))then
193     ! ori tvp(i,k)=tvp(i,k)*(1.0-qnk(i)+ep(i,k)*clw(i,k))
194     ! oric print*,'i,k,tvp(i,k),qnk(i),ep(i,k),clw(i,k)'
195     ! oric print*, i,k,tvp(i,k),qnk(i),ep(i,k),clw(i,k)
196     ! ori endif
197     ! ori 330 continue
198     ! ori 340 continue
199    
200     ! ori do 350 i=1,ncum
201     ! ori tvp(i,nlp)=tvp(i,nl)-(gz(i,nlp)-gz(i,nl))/cpd
202     ! ori 350 continue
203    
204     do 350 i=1,ncum ! convect3
205     tp(i,nlp)=tp(i,nl) ! convect3
206     350 continue ! convect3
207     !
208     !=====================================================================
209     ! --- EFFECTIVE VERTICAL PROFILE OF BUOYANCY (convect3 only):
210     !=====================================================================
211    
212     !-- this is for convect3 only:
213    
214     ! first estimate of buoyancy:
215    
216     do 500 i=1,ncum
217     do 501 k=1,nl
218     buoy(i,k)=tvp(i,k)-tv(i,k)
219     501 continue
220     500 continue
221    
222     ! set buoyancy=buoybase for all levels below base
223     ! for safety, set buoy(icb)=buoybase
224    
225     do 505 i=1,ncum
226     do 506 k=1,nl
227     if((k.ge.icb(i)).and.(k.le.nl).and.(p(i,k).ge.pbase(i)))then
228     buoy(i,k)=buoybase(i)
229     endif
230     506 continue
231     buoy(icb(i),k)=buoybase(i)
232     505 continue
233    
234     !-- end convect3
235    
236     !=====================================================================
237     ! --- FIND THE FIRST MODEL LEVEL (INB) ABOVE THE PARCEL'S
238     ! --- LEVEL OF NEUTRAL BUOYANCY
239     !=====================================================================
240     !
241     !-- this is for convect3 only:
242    
243     do 510 i=1,ncum
244     inb(i)=nl-1
245     510 continue
246    
247     do 530 i=1,ncum
248     do 535 k=1,nl-1
249     if ((k.ge.icb(i)).and.(buoy(i,k).lt.dtovsh)) then
250     inb(i)=MIN(inb(i),k)
251     endif
252     535 continue
253     530 continue
254    
255     !-- end convect3
256    
257     ! ori do 510 i=1,ncum
258     ! ori cape(i)=0.0
259     ! ori capem(i)=0.0
260     ! ori inb(i)=icb(i)+1
261     ! ori inb1(i)=inb(i)
262     ! ori 510 continue
263     !
264     ! Originial Code
265     !
266     ! do 530 k=minorig+1,nl-1
267     ! do 520 i=1,ncum
268     ! if(k.ge.(icb(i)+1))then
269     ! by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k)
270     ! byp=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1)
271     ! cape(i)=cape(i)+by
272     ! if(by.ge.0.0)inb1(i)=k+1
273     ! if(cape(i).gt.0.0)then
274     ! inb(i)=k+1
275     ! capem(i)=cape(i)
276     ! endif
277     ! endif
278     !520 continue
279     !530 continue
280     ! do 540 i=1,ncum
281     ! byp=(tvp(i,nl)-tv(i,nl))*dph(i,nl)/p(i,nl)
282     ! cape(i)=capem(i)+byp
283     ! defrac=capem(i)-cape(i)
284     ! defrac=max(defrac,0.001)
285     ! frac(i)=-cape(i)/defrac
286     ! frac(i)=min(frac(i),1.0)
287     ! frac(i)=max(frac(i),0.0)
288     !540 continue
289     !
290     ! K Emanuel fix
291     !
292     ! call zilch(byp,ncum)
293     ! do 530 k=minorig+1,nl-1
294     ! do 520 i=1,ncum
295     ! if(k.ge.(icb(i)+1))then
296     ! by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k)
297     ! cape(i)=cape(i)+by
298     ! if(by.ge.0.0)inb1(i)=k+1
299     ! if(cape(i).gt.0.0)then
300     ! inb(i)=k+1
301     ! capem(i)=cape(i)
302     ! byp(i)=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1)
303     ! endif
304     ! endif
305     !520 continue
306     !530 continue
307     ! do 540 i=1,ncum
308     ! inb(i)=max(inb(i),inb1(i))
309     ! cape(i)=capem(i)+byp(i)
310     ! defrac=capem(i)-cape(i)
311     ! defrac=max(defrac,0.001)
312     ! frac(i)=-cape(i)/defrac
313     ! frac(i)=min(frac(i),1.0)
314     ! frac(i)=max(frac(i),0.0)
315     !540 continue
316     !
317     ! J Teixeira fix
318     !
319     ! ori call zilch(byp,ncum)
320     ! ori do 515 i=1,ncum
321     ! ori lcape(i)=.true.
322     ! ori 515 continue
323     ! ori do 530 k=minorig+1,nl-1
324     ! ori do 520 i=1,ncum
325     ! ori if(cape(i).lt.0.0)lcape(i)=.false.
326     ! ori if((k.ge.(icb(i)+1)).and.lcape(i))then
327     ! ori by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k)
328     ! ori byp(i)=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1)
329     ! ori cape(i)=cape(i)+by
330     ! ori if(by.ge.0.0)inb1(i)=k+1
331     ! ori if(cape(i).gt.0.0)then
332     ! ori inb(i)=k+1
333     ! ori capem(i)=cape(i)
334     ! ori endif
335     ! ori endif
336     ! ori 520 continue
337     ! ori 530 continue
338     ! ori do 540 i=1,ncum
339     ! ori cape(i)=capem(i)+byp(i)
340     ! ori defrac=capem(i)-cape(i)
341     ! ori defrac=max(defrac,0.001)
342     ! ori frac(i)=-cape(i)/defrac
343     ! ori frac(i)=min(frac(i),1.0)
344     ! ori frac(i)=max(frac(i),0.0)
345     ! ori 540 continue
346     !
347     !=====================================================================
348     ! --- CALCULATE LIQUID WATER STATIC ENERGY OF LIFTED PARCEL
349     !=====================================================================
350     !
351     !ym do i=1,ncum*nlp
352     !ym hp(i,1)=h(i,1)
353     !ym enddo
354    
355     do k=1,nlp
356     do i=1,ncum
357     hp(i,k)=h(i,k)
358     enddo
359     enddo
360    
361     do 600 k=minorig+1,nl
362     do 590 i=1,ncum
363     if((k.ge.icb(i)).and.(k.le.inb(i)))then
364     hp(i,k)=h(i,nk(i))+(lv(i,k)+(cpd-cpv)*t(i,k))*ep(i,k)*clw(i,k)
365     endif
366     590 continue
367     600 continue
368    
369     return
370     end

  ViewVC Help
Powered by ViewVC 1.1.21