/[lmdze]/trunk/Sources/phylmd/CV3_routines/cv3_unsat.f
ViewVC logotype

Annotation of /trunk/Sources/phylmd/CV3_routines/cv3_unsat.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 47 - (hide annotations)
Fri Jul 1 15:00:48 2011 UTC (12 years, 11 months ago) by guez
Original Path: trunk/libf/phylmd/CV3_routines/cv3_unsat.f90
File size: 10213 byte(s)
Split "thermcell.f" and "cv3_routines.f".
Removed copies of files that are now in "L_util".
Moved "mva9" and "diagetpq" to their own files.
Unified variable names across procedures.

1 guez 47
2    
3     SUBROUTINE cv3_unsat(nloc,ncum,nd,na,ntra,icb,inb &
4     ,t,rr,rs,gz,u,v,tra,p,ph &
5     ,th,tv,lv,cpn,ep,sigp,clw &
6     ,m,ment,elij,delt,plcl &
7     ,mp,rp,up,vp,trap,wt,water,evap,b)
8     use cvparam3
9     use cvthermo
10     use cvflag
11     implicit none
12    
13    
14    
15     ! inputs:
16     integer ncum, nd, na, ntra, nloc
17     integer icb(nloc), inb(nloc)
18     real, intent(in):: delt
19     real plcl(nloc)
20     real t(nloc,nd), rr(nloc,nd), rs(nloc,nd)
21     real u(nloc,nd), v(nloc,nd)
22     real tra(nloc,nd,ntra)
23     real p(nloc,nd), ph(nloc,nd+1)
24     real th(nloc,na), gz(nloc,na)
25     real lv(nloc,na), ep(nloc,na), sigp(nloc,na), clw(nloc,na)
26     real cpn(nloc,na), tv(nloc,na)
27     real m(nloc,na), ment(nloc,na,na), elij(nloc,na,na)
28    
29     ! outputs:
30     real mp(nloc,na), rp(nloc,na), up(nloc,na), vp(nloc,na)
31     real water(nloc,na), evap(nloc,na), wt(nloc,na)
32     real trap(nloc,na,ntra)
33     real b(nloc,na)
34    
35     ! local variables
36     integer i,j,k,il,num1
37     real tinv, delti
38     real awat, afac, afac1, afac2, bfac
39     real pr1, pr2, sigt, b6, c6, revap, tevap, delth
40     real amfac, amp2, xf, tf, fac2, ur, sru, fac, d, af, bf
41     real ampmax
42     real lvcp(nloc,na)
43     real wdtrain(nloc)
44     logical lwork(nloc)
45    
46    
47     !------------------------------------------------------
48    
49     delti = 1./delt
50     tinv=1./3.
51    
52     mp(:,:)=0.
53    
54     do i=1,nl
55     do il=1,ncum
56     mp(il,i)=0.0
57     rp(il,i)=rr(il,i)
58     up(il,i)=u(il,i)
59     vp(il,i)=v(il,i)
60     wt(il,i)=0.001
61     water(il,i)=0.0
62     evap(il,i)=0.0
63     b(il,i)=0.0
64     lvcp(il,i)=lv(il,i)/cpn(il,i)
65     enddo
66     enddo
67    
68     ! do k=1,ntra
69     ! do i=1,nd
70     ! do il=1,ncum
71     ! trap(il,i,k)=tra(il,i,k)
72     ! enddo
73     ! enddo
74     ! enddo
75    
76     !
77     ! *** check whether ep(inb)=0, if so, skip precipitating ***
78     ! *** downdraft calculation ***
79     !
80    
81     do il=1,ncum
82     lwork(il)=.TRUE.
83     if(ep(il,inb(il)).lt.0.0001)lwork(il)=.FALSE.
84     enddo
85    
86     call zilch(wdtrain,ncum)
87    
88     DO 400 i=nl+1,1,-1
89    
90     num1=0
91     do il=1,ncum
92     if ( i.le.inb(il) .and. lwork(il) ) num1=num1+1
93     enddo
94     if (num1.le.0) goto 400
95    
96     !
97     ! *** integrate liquid water equation to find condensed water ***
98     ! *** and condensed water flux ***
99     !
100    
101     !
102     ! *** begin downdraft loop ***
103     !
104    
105     !
106     ! *** calculate detrained precipitation ***
107     !
108     do il=1,ncum
109     if (i.le.inb(il) .and. lwork(il)) then
110     if (cvflag_grav) then
111     wdtrain(il)=grav*ep(il,i)*m(il,i)*clw(il,i)
112     else
113     wdtrain(il)=10.0*ep(il,i)*m(il,i)*clw(il,i)
114     endif
115     endif
116     enddo
117    
118     if(i.gt.1)then
119     do 320 j=1,i-1
120     do il=1,ncum
121     if (i.le.inb(il) .and. lwork(il)) then
122     awat=elij(il,j,i)-(1.-ep(il,i))*clw(il,i)
123     awat=amax1(awat,0.0)
124     if (cvflag_grav) then
125     wdtrain(il)=wdtrain(il)+grav*awat*ment(il,j,i)
126     else
127     wdtrain(il)=wdtrain(il)+10.0*awat*ment(il,j,i)
128     endif
129     endif
130     enddo
131     320 continue
132     endif
133    
134     !
135     ! *** find rain water and evaporation using provisional ***
136     ! *** estimates of rp(i)and rp(i-1) ***
137     !
138    
139     do 999 il=1,ncum
140    
141     if (i.le.inb(il) .and. lwork(il)) then
142    
143     wt(il,i)=45.0
144    
145     if(i.lt.inb(il))then
146     rp(il,i)=rp(il,i+1) &
147     +(cpd*(t(il,i+1)-t(il,i))+gz(il,i+1)-gz(il,i))/lv(il,i)
148     rp(il,i)=0.5*(rp(il,i)+rr(il,i))
149     endif
150     rp(il,i)=amax1(rp(il,i),0.0)
151     rp(il,i)=amin1(rp(il,i),rs(il,i))
152     rp(il,inb(il))=rr(il,inb(il))
153    
154     if(i.eq.1)then
155     afac=p(il,1)*(rs(il,1)-rp(il,1))/(1.0e4+2000.0*p(il,1)*rs(il,1))
156     else
157     rp(il,i-1)=rp(il,i) &
158     +(cpd*(t(il,i)-t(il,i-1))+gz(il,i)-gz(il,i-1))/lv(il,i)
159     rp(il,i-1)=0.5*(rp(il,i-1)+rr(il,i-1))
160     rp(il,i-1)=amin1(rp(il,i-1),rs(il,i-1))
161     rp(il,i-1)=amax1(rp(il,i-1),0.0)
162     afac1=p(il,i)*(rs(il,i)-rp(il,i))/(1.0e4+2000.0*p(il,i)*rs(il,i))
163     afac2=p(il,i-1)*(rs(il,i-1)-rp(il,i-1)) &
164     /(1.0e4+2000.0*p(il,i-1)*rs(il,i-1))
165     afac=0.5*(afac1+afac2)
166     endif
167     if(i.eq.inb(il))afac=0.0
168     afac=amax1(afac,0.0)
169     bfac=1./(sigd*wt(il,i))
170     !
171     !jyg1
172     !cc sigt=1.0
173     !cc if(i.ge.icb)sigt=sigp(i)
174     ! prise en compte de la variation progressive de sigt dans
175     ! les couches icb et icb-1:
176     ! pour plcl<ph(i+1), pr1=0 & pr2=1
177     ! pour plcl>ph(i), pr1=1 & pr2=0
178     ! pour ph(i+1)<plcl<ph(i), pr1 est la proportion a cheval
179     ! sur le nuage, et pr2 est la proportion sous la base du
180     ! nuage.
181     pr1=(plcl(il)-ph(il,i+1))/(ph(il,i)-ph(il,i+1))
182     pr1=max(0.,min(1.,pr1))
183     pr2=(ph(il,i)-plcl(il))/(ph(il,i)-ph(il,i+1))
184     pr2=max(0.,min(1.,pr2))
185     sigt=sigp(il,i)*pr1+pr2
186     !jyg2
187     !
188     b6=bfac*50.*sigd*(ph(il,i)-ph(il,i+1))*sigt*afac
189     c6=water(il,i+1)+bfac*wdtrain(il) &
190     -50.*sigd*bfac*(ph(il,i)-ph(il,i+1))*evap(il,i+1)
191     if(c6.gt.0.0)then
192     revap=0.5*(-b6+sqrt(b6*b6+4.*c6))
193     evap(il,i)=sigt*afac*revap
194     water(il,i)=revap*revap
195     else
196     evap(il,i)=-evap(il,i+1) &
197     +0.02*(wdtrain(il)+sigd*wt(il,i)*water(il,i+1)) &
198     /(sigd*(ph(il,i)-ph(il,i+1)))
199     end if
200     !
201     ! *** calculate precipitating downdraft mass flux under ***
202     ! *** hydrostatic approximation ***
203     !
204     if (i.ne.1) then
205    
206     tevap=amax1(0.0,evap(il,i))
207     delth=amax1(0.001,(th(il,i)-th(il,i-1)))
208     if (cvflag_grav) then
209     mp(il,i)=100.*ginv*lvcp(il,i)*sigd*tevap &
210     *(p(il,i-1)-p(il,i))/delth
211     else
212     mp(il,i)=10.*lvcp(il,i)*sigd*tevap*(p(il,i-1)-p(il,i))/delth
213     endif
214     !
215     ! *** if hydrostatic assumption fails, ***
216     ! *** solve cubic difference equation for downdraft theta ***
217     ! *** and mass flux from two simultaneous differential eqns ***
218     !
219     amfac=sigd*sigd*70.0*ph(il,i)*(p(il,i-1)-p(il,i)) &
220     *(th(il,i)-th(il,i-1))/(tv(il,i)*th(il,i))
221     amp2=abs(mp(il,i+1)*mp(il,i+1)-mp(il,i)*mp(il,i))
222     if(amp2.gt.(0.1*amfac))then
223     xf=100.0*sigd*sigd*sigd*(ph(il,i)-ph(il,i+1))
224     tf=b(il,i)-5.0*(th(il,i)-th(il,i-1))*t(il,i) &
225     /(lvcp(il,i)*sigd*th(il,i))
226     af=xf*tf+mp(il,i+1)*mp(il,i+1)*tinv
227     bf=2.*(tinv*mp(il,i+1))**3+tinv*mp(il,i+1)*xf*tf &
228     +50.*(p(il,i-1)-p(il,i))*xf*tevap
229     fac2=1.0
230     if(bf.lt.0.0)fac2=-1.0
231     bf=abs(bf)
232     ur=0.25*bf*bf-af*af*af*tinv*tinv*tinv
233     if(ur.ge.0.0)then
234     sru=sqrt(ur)
235     fac=1.0
236     if((0.5*bf-sru).lt.0.0)fac=-1.0
237     mp(il,i)=mp(il,i+1)*tinv+(0.5*bf+sru)**tinv &
238     +fac*(abs(0.5*bf-sru))**tinv
239     else
240     d=atan(2.*sqrt(-ur)/(bf+1.0e-28))
241     if(fac2.lt.0.0)d=3.14159-d
242     mp(il,i)=mp(il,i+1)*tinv+2.*sqrt(af*tinv)*cos(d*tinv)
243     endif
244     mp(il,i)=amax1(0.0,mp(il,i))
245    
246     if (cvflag_grav) then
247     !jyg : il y a vraisemblablement une erreur dans la ligne 2 suivante:
248     ! il faut diviser par (mp(il,i)*sigd*grav) et non par (mp(il,i)+sigd*0.1).
249     ! Et il faut bien revoir les facteurs 100.
250     b(il,i-1)=b(il,i)+100.0*(p(il,i-1)-p(il,i))*tevap &
251     /(mp(il,i)+sigd*0.1) &
252     -10.0*(th(il,i)-th(il,i-1))*t(il,i)/(lvcp(il,i)*sigd*th(il,i))
253     else
254     b(il,i-1)=b(il,i)+100.0*(p(il,i-1)-p(il,i))*tevap &
255     /(mp(il,i)+sigd*0.1) &
256     -10.0*(th(il,i)-th(il,i-1))*t(il,i)/(lvcp(il,i)*sigd*th(il,i))
257     endif
258     b(il,i-1)=amax1(b(il,i-1),0.0)
259     endif
260     !
261     ! *** limit magnitude of mp(i) to meet cfl condition ***
262     !
263     ampmax=2.0*(ph(il,i)-ph(il,i+1))*delti
264     amp2=2.0*(ph(il,i-1)-ph(il,i))*delti
265     ampmax=amin1(ampmax,amp2)
266     mp(il,i)=amin1(mp(il,i),ampmax)
267     !
268     ! *** force mp to decrease linearly to zero ***
269     ! *** between cloud base and the surface ***
270     !
271     if(p(il,i).gt.p(il,icb(il)))then
272     mp(il,i)=mp(il,icb(il))*(p(il,1)-p(il,i))/(p(il,1)-p(il,icb(il)))
273     endif
274    
275     360 continue
276     endif ! i.eq.1
277     !
278     ! *** find mixing ratio of precipitating downdraft ***
279     !
280    
281     if (i.ne.inb(il)) then
282    
283     rp(il,i)=rr(il,i)
284    
285     if(mp(il,i).gt.mp(il,i+1))then
286    
287     if (cvflag_grav) then
288     rp(il,i)=rp(il,i+1)*mp(il,i+1)+rr(il,i)*(mp(il,i)-mp(il,i+1)) &
289     +100.*ginv*0.5*sigd*(ph(il,i)-ph(il,i+1)) &
290     *(evap(il,i+1)+evap(il,i))
291     else
292     rp(il,i)=rp(il,i+1)*mp(il,i+1)+rr(il,i)*(mp(il,i)-mp(il,i+1)) &
293     +5.*sigd*(ph(il,i)-ph(il,i+1)) &
294     *(evap(il,i+1)+evap(il,i))
295     endif
296     rp(il,i)=rp(il,i)/mp(il,i)
297     up(il,i)=up(il,i+1)*mp(il,i+1)+u(il,i)*(mp(il,i)-mp(il,i+1))
298     up(il,i)=up(il,i)/mp(il,i)
299     vp(il,i)=vp(il,i+1)*mp(il,i+1)+v(il,i)*(mp(il,i)-mp(il,i+1))
300     vp(il,i)=vp(il,i)/mp(il,i)
301    
302     ! do j=1,ntra
303     ! trap(il,i,j)=trap(il,i+1,j)*mp(il,i+1)
304     !testmaf : +trap(il,i,j)*(mp(il,i)-mp(il,i+1))
305     ! : +tra(il,i,j)*(mp(il,i)-mp(il,i+1))
306     ! trap(il,i,j)=trap(il,i,j)/mp(il,i)
307     ! end do
308    
309     else
310    
311     if(mp(il,i+1).gt.1.0e-16)then
312     if (cvflag_grav) then
313     rp(il,i)=rp(il,i+1) &
314     +100.*ginv*0.5*sigd*(ph(il,i)-ph(il,i+1)) &
315     *(evap(il,i+1)+evap(il,i))/mp(il,i+1)
316     else
317     rp(il,i)=rp(il,i+1) &
318     +5.*sigd*(ph(il,i)-ph(il,i+1)) &
319     *(evap(il,i+1)+evap(il,i))/mp(il,i+1)
320     endif
321     up(il,i)=up(il,i+1)
322     vp(il,i)=vp(il,i+1)
323    
324     ! do j=1,ntra
325     ! trap(il,i,j)=trap(il,i+1,j)
326     ! end do
327    
328     endif
329     endif
330     rp(il,i)=amin1(rp(il,i),rs(il,i))
331     rp(il,i)=amax1(rp(il,i),0.0)
332    
333     endif
334     endif
335     999 continue
336    
337     400 continue
338    
339     return
340     end

  ViewVC Help
Powered by ViewVC 1.1.21