/[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 178 - (hide annotations)
Fri Mar 11 18:47:26 2016 UTC (8 years, 2 months ago) by guez
Original Path: trunk/Sources/phylmd/CV3_routines/cv3_undilute2.f
File size: 11364 byte(s)
Moved variables date0, deltat, datasz_max, ncvar_ids, point, buff_pos,
buffer, regular from module histcom_var to modules where they are
defined.

Removed procedure ioipslmpp, useless for a sequential program.

Added argument datasz_max to histwrite_real (to avoid circular
dependency with histwrite).

Removed useless variables and computations everywhere.

Changed real litteral constants from default kind to double precision
in lwb, lwu, lwvn, sw1s, swtt, swtt1, swu.

Removed unused arguments: paer of sw, sw1s, sw2s, swclr; pcldsw of
sw1s, sw2s; pdsig, prayl of swr; co2_ppm of clmain, clqh; tsol of
transp_lay; nsrf of screenp; kcrit and kknu of gwstress; pstd of
orosetup.

Added output of relative humidity.

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 178 real tg,qg,ahg,alv,s,tc,es,denom
48 guez 150 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     ! convect3: (qg utilise au lieu du vrai mixing ratio rg):
145     tvp(i,k)=tp(i,k)*(1.+qg/eps-qnk(i)) ! whole thing
146     endif
147     290 continue
148     300 continue
149     !
150     !=====================================================================
151     ! --- SET THE PRECIPITATION EFFICIENCIES AND THE FRACTION OF
152     ! --- PRECIPITATION FALLING OUTSIDE OF CLOUD
153     ! --- THESE MAY BE FUNCTIONS OF TP(I), P(I) AND CLW(I)
154     !=====================================================================
155     !
156     ! ori do 320 k=minorig+1,nl
157     do 320 k=1,nl ! convect3
158     do 310 i=1,ncum
159     pden=ptcrit-pbcrit
160     ep(i,k)=(plcl(i)-p(i,k)-pbcrit)/pden*epmax
161     ep(i,k)=amax1(ep(i,k),0.0)
162     ep(i,k)=amin1(ep(i,k),epmax)
163     sigp(i,k)=spfac
164     ! ori if(k.ge.(nk(i)+1))then
165     ! ori tca=tp(i,k)-t0
166     ! ori if(tca.ge.0.0)then
167     ! ori elacrit=elcrit
168     ! ori else
169     ! ori elacrit=elcrit*(1.0-tca/tlcrit)
170     ! ori endif
171     ! ori elacrit=max(elacrit,0.0)
172     ! ori ep(i,k)=1.0-elacrit/max(clw(i,k),1.0e-8)
173     ! ori ep(i,k)=max(ep(i,k),0.0 )
174     ! ori ep(i,k)=min(ep(i,k),1.0 )
175     ! ori sigp(i,k)=sigs
176     ! ori endif
177     310 continue
178     320 continue
179     !
180     !=====================================================================
181     ! --- CALCULATE VIRTUAL TEMPERATURE AND LIFTED PARCEL
182     ! --- VIRTUAL TEMPERATURE
183     !=====================================================================
184     !
185     ! dans convect3, tvp est calcule en une seule fois, et sans retirer
186     ! l'eau condensee (~> reversible CAPE)
187     !
188     ! ori do 340 k=minorig+1,nl
189     ! ori do 330 i=1,ncum
190     ! ori if(k.ge.(icb(i)+1))then
191     ! ori tvp(i,k)=tvp(i,k)*(1.0-qnk(i)+ep(i,k)*clw(i,k))
192     ! oric print*,'i,k,tvp(i,k),qnk(i),ep(i,k),clw(i,k)'
193     ! oric print*, i,k,tvp(i,k),qnk(i),ep(i,k),clw(i,k)
194     ! ori endif
195     ! ori 330 continue
196     ! ori 340 continue
197    
198     ! ori do 350 i=1,ncum
199     ! ori tvp(i,nlp)=tvp(i,nl)-(gz(i,nlp)-gz(i,nl))/cpd
200     ! ori 350 continue
201    
202     do 350 i=1,ncum ! convect3
203     tp(i,nlp)=tp(i,nl) ! convect3
204     350 continue ! convect3
205     !
206     !=====================================================================
207     ! --- EFFECTIVE VERTICAL PROFILE OF BUOYANCY (convect3 only):
208     !=====================================================================
209    
210     !-- this is for convect3 only:
211    
212     ! first estimate of buoyancy:
213    
214     do 500 i=1,ncum
215     do 501 k=1,nl
216     buoy(i,k)=tvp(i,k)-tv(i,k)
217     501 continue
218     500 continue
219    
220     ! set buoyancy=buoybase for all levels below base
221     ! for safety, set buoy(icb)=buoybase
222    
223     do 505 i=1,ncum
224     do 506 k=1,nl
225     if((k.ge.icb(i)).and.(k.le.nl).and.(p(i,k).ge.pbase(i)))then
226     buoy(i,k)=buoybase(i)
227     endif
228     506 continue
229     buoy(icb(i),k)=buoybase(i)
230     505 continue
231    
232     !-- end convect3
233    
234     !=====================================================================
235     ! --- FIND THE FIRST MODEL LEVEL (INB) ABOVE THE PARCEL'S
236     ! --- LEVEL OF NEUTRAL BUOYANCY
237     !=====================================================================
238     !
239     !-- this is for convect3 only:
240    
241     do 510 i=1,ncum
242     inb(i)=nl-1
243     510 continue
244    
245     do 530 i=1,ncum
246     do 535 k=1,nl-1
247     if ((k.ge.icb(i)).and.(buoy(i,k).lt.dtovsh)) then
248     inb(i)=MIN(inb(i),k)
249     endif
250     535 continue
251     530 continue
252    
253     !-- end convect3
254    
255     ! ori do 510 i=1,ncum
256     ! ori cape(i)=0.0
257     ! ori capem(i)=0.0
258     ! ori inb(i)=icb(i)+1
259     ! ori inb1(i)=inb(i)
260     ! ori 510 continue
261     !
262     ! Originial Code
263     !
264     ! do 530 k=minorig+1,nl-1
265     ! do 520 i=1,ncum
266     ! if(k.ge.(icb(i)+1))then
267     ! by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k)
268     ! byp=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1)
269     ! cape(i)=cape(i)+by
270     ! if(by.ge.0.0)inb1(i)=k+1
271     ! if(cape(i).gt.0.0)then
272     ! inb(i)=k+1
273     ! capem(i)=cape(i)
274     ! endif
275     ! endif
276     !520 continue
277     !530 continue
278     ! do 540 i=1,ncum
279     ! byp=(tvp(i,nl)-tv(i,nl))*dph(i,nl)/p(i,nl)
280     ! cape(i)=capem(i)+byp
281     ! defrac=capem(i)-cape(i)
282     ! defrac=max(defrac,0.001)
283     ! frac(i)=-cape(i)/defrac
284     ! frac(i)=min(frac(i),1.0)
285     ! frac(i)=max(frac(i),0.0)
286     !540 continue
287     !
288     ! K Emanuel fix
289     !
290     ! call zilch(byp,ncum)
291     ! do 530 k=minorig+1,nl-1
292     ! do 520 i=1,ncum
293     ! if(k.ge.(icb(i)+1))then
294     ! by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k)
295     ! cape(i)=cape(i)+by
296     ! if(by.ge.0.0)inb1(i)=k+1
297     ! if(cape(i).gt.0.0)then
298     ! inb(i)=k+1
299     ! capem(i)=cape(i)
300     ! byp(i)=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1)
301     ! endif
302     ! endif
303     !520 continue
304     !530 continue
305     ! do 540 i=1,ncum
306     ! inb(i)=max(inb(i),inb1(i))
307     ! cape(i)=capem(i)+byp(i)
308     ! defrac=capem(i)-cape(i)
309     ! defrac=max(defrac,0.001)
310     ! frac(i)=-cape(i)/defrac
311     ! frac(i)=min(frac(i),1.0)
312     ! frac(i)=max(frac(i),0.0)
313     !540 continue
314     !
315     ! J Teixeira fix
316     !
317     ! ori call zilch(byp,ncum)
318     ! ori do 515 i=1,ncum
319     ! ori lcape(i)=.true.
320     ! ori 515 continue
321     ! ori do 530 k=minorig+1,nl-1
322     ! ori do 520 i=1,ncum
323     ! ori if(cape(i).lt.0.0)lcape(i)=.false.
324     ! ori if((k.ge.(icb(i)+1)).and.lcape(i))then
325     ! ori by=(tvp(i,k)-tv(i,k))*dph(i,k)/p(i,k)
326     ! ori byp(i)=(tvp(i,k+1)-tv(i,k+1))*dph(i,k+1)/p(i,k+1)
327     ! ori cape(i)=cape(i)+by
328     ! ori if(by.ge.0.0)inb1(i)=k+1
329     ! ori if(cape(i).gt.0.0)then
330     ! ori inb(i)=k+1
331     ! ori capem(i)=cape(i)
332     ! ori endif
333     ! ori endif
334     ! ori 520 continue
335     ! ori 530 continue
336     ! ori do 540 i=1,ncum
337     ! ori cape(i)=capem(i)+byp(i)
338     ! ori defrac=capem(i)-cape(i)
339     ! ori defrac=max(defrac,0.001)
340     ! ori frac(i)=-cape(i)/defrac
341     ! ori frac(i)=min(frac(i),1.0)
342     ! ori frac(i)=max(frac(i),0.0)
343     ! ori 540 continue
344     !
345     !=====================================================================
346     ! --- CALCULATE LIQUID WATER STATIC ENERGY OF LIFTED PARCEL
347     !=====================================================================
348     !
349     !ym do i=1,ncum*nlp
350     !ym hp(i,1)=h(i,1)
351     !ym enddo
352    
353     do k=1,nlp
354     do i=1,ncum
355     hp(i,k)=h(i,k)
356     enddo
357     enddo
358    
359     do 600 k=minorig+1,nl
360     do 590 i=1,ncum
361     if((k.ge.icb(i)).and.(k.le.inb(i)))then
362     hp(i,k)=h(i,nk(i))+(lv(i,k)+(cpd-cpv)*t(i,k))*ep(i,k)*clw(i,k)
363     endif
364     590 continue
365     600 continue
366    
367     return
368     end

  ViewVC Help
Powered by ViewVC 1.1.21