/[lmdze]/trunk/phylmd/CV_routines/cv_unsat.f
ViewVC logotype

Annotation of /trunk/phylmd/CV_routines/cv_unsat.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 52 - (hide annotations)
Fri Sep 23 12:28:01 2011 UTC (12 years, 8 months ago) by guez
Original Path: trunk/libf/phylmd/CV_routines/cv_unsat.f90
File size: 6185 byte(s)
Split "conflx.f" into single-procedure files in directory "Conflx".

Split "cv_routines.f" into single-procedure files in directory
"CV_routines". Made module "cvparam" from included file
"cvparam.h". No included file other than "netcdf.inc" left in LMDZE.

1 guez 52
2     SUBROUTINE cv_unsat(nloc,ncum,nd,inb,t,q,qs,gz,u,v,p,ph &
3     ,h,lv,ep,sigp,clw,m,ment,elij &
4     ,iflag,mp,qp,up,vp,wt,water,evap)
5     use cvthermo
6     use cvparam
7     implicit none
8    
9    
10    
11     ! inputs:
12     integer ncum, nd, nloc
13     integer inb(nloc)
14     real t(nloc,nd), q(nloc,nd), qs(nloc,nd)
15     real gz(nloc,nd), u(nloc,nd), v(nloc,nd)
16     real p(nloc,nd), ph(nloc,nd+1), h(nloc,nd)
17     real lv(nloc,nd), ep(nloc,nd), sigp(nloc,nd), clw(nloc,nd)
18     real m(nloc,nd), ment(nloc,nd,nd), elij(nloc,nd,nd)
19    
20     ! outputs:
21     integer iflag(nloc) ! also an input
22     real mp(nloc,nd), qp(nloc,nd), up(nloc,nd), vp(nloc,nd)
23     real water(nloc,nd), evap(nloc,nd), wt(nloc,nd)
24    
25     ! local variables:
26     integer i,j,k,ij,num1
27     integer jtt(nloc)
28     real awat, coeff, qsm, afac, sigt, b6, c6, revap
29     real dhdp, fac, qstm, rat
30     real wdtrain(nloc)
31     logical lwork(nloc)
32    
33     !=====================================================================
34     ! --- PRECIPITATING DOWNDRAFT CALCULATION
35     !=====================================================================
36     !
37     ! Initializations:
38     !
39     do i = 1, ncum
40     do k = 1, nl+1
41     wt(i,k) = omtsnow
42     mp(i,k) = 0.0
43     evap(i,k) = 0.0
44     water(i,k) = 0.0
45     enddo
46     enddo
47    
48     do 420 i=1,ncum
49     qp(i,1)=q(i,1)
50     up(i,1)=u(i,1)
51     vp(i,1)=v(i,1)
52     420 continue
53    
54     do 440 k=2,nl+1
55     do 430 i=1,ncum
56     qp(i,k)=q(i,k-1)
57     up(i,k)=u(i,k-1)
58     vp(i,k)=v(i,k-1)
59     430 continue
60     440 continue
61    
62    
63     ! *** Check whether ep(inb)=0, if so, skip precipitating ***
64     ! *** downdraft calculation ***
65     !
66     !
67     ! *** Integrate liquid water equation to find condensed water ***
68     ! *** and condensed water flux ***
69     !
70     !
71     do 890 i=1,ncum
72     jtt(i)=2
73     if(ep(i,inb(i)).le.0.0001)iflag(i)=2
74     if(iflag(i).eq.0)then
75     lwork(i)=.true.
76     else
77     lwork(i)=.false.
78     endif
79     890 continue
80     !
81     ! *** Begin downdraft loop ***
82     !
83     !
84     call zilch(wdtrain,ncum)
85     do 899 i=nl+1,1,-1
86     !
87     num1=0
88     do 879 ij=1,ncum
89     if((i.le.inb(ij)).and.lwork(ij))num1=num1+1
90     879 continue
91     if(num1.le.0)go to 899
92     !
93     !
94     ! *** Calculate detrained precipitation ***
95     !
96     do 891 ij=1,ncum
97     if((i.le.inb(ij)).and.(lwork(ij)))then
98     wdtrain(ij)=g*ep(ij,i)*m(ij,i)*clw(ij,i)
99     endif
100     891 continue
101     !
102     if(i.gt.1)then
103     do 893 j=1,i-1
104     do 892 ij=1,ncum
105     if((i.le.inb(ij)).and.(lwork(ij)))then
106     awat=elij(ij,j,i)-(1.-ep(ij,i))*clw(ij,i)
107     awat=max(0.0,awat)
108     wdtrain(ij)=wdtrain(ij)+g*awat*ment(ij,j,i)
109     endif
110     892 continue
111     893 continue
112     endif
113     !
114     ! *** Find rain water and evaporation using provisional ***
115     ! *** estimates of qp(i)and qp(i-1) ***
116     !
117     !
118     ! *** Value of terminal velocity and coeffecient of evaporation for snow ***
119     !
120     do 894 ij=1,ncum
121     if((i.le.inb(ij)).and.(lwork(ij)))then
122     coeff=coeffs
123     wt(ij,i)=omtsnow
124     !
125     ! *** Value of terminal velocity and coeffecient of evaporation for rain ***
126     !
127     if(t(ij,i).gt.273.0)then
128     coeff=coeffr
129     wt(ij,i)=omtrain
130     endif
131     qsm=0.5*(q(ij,i)+qp(ij,i+1))
132     afac=coeff*ph(ij,i)*(qs(ij,i)-qsm) &
133     /(1.0e4+2.0e3*ph(ij,i)*qs(ij,i))
134     afac=max(afac,0.0)
135     sigt=sigp(ij,i)
136     sigt=max(0.0,sigt)
137     sigt=min(1.0,sigt)
138     b6=100.*(ph(ij,i)-ph(ij,i+1))*sigt*afac/wt(ij,i)
139     c6=(water(ij,i+1)*wt(ij,i+1)+wdtrain(ij)/sigd)/wt(ij,i)
140     revap=0.5*(-b6+sqrt(b6*b6+4.*c6))
141     evap(ij,i)=sigt*afac*revap
142     water(ij,i)=revap*revap
143     !
144     ! *** Calculate precipitating downdraft mass flux under ***
145     ! *** hydrostatic approximation ***
146     !
147     if(i.gt.1)then
148     dhdp=(h(ij,i)-h(ij,i-1))/(p(ij,i-1)-p(ij,i))
149     dhdp=max(dhdp,10.0)
150     mp(ij,i)=100.*ginv*lv(ij,i)*sigd*evap(ij,i)/dhdp
151     mp(ij,i)=max(mp(ij,i),0.0)
152     !
153     ! *** Add small amount of inertia to downdraft ***
154     !
155     fac=20.0/(ph(ij,i-1)-ph(ij,i))
156     mp(ij,i)=(fac*mp(ij,i+1)+mp(ij,i))/(1.+fac)
157     !
158     ! *** Force mp to decrease linearly to zero ***
159     ! *** between about 950 mb and the surface ***
160     !
161     if(p(ij,i).gt.(0.949*p(ij,1)))then
162     jtt(ij)=max(jtt(ij),i)
163     mp(ij,i)=mp(ij,jtt(ij))*(p(ij,1)-p(ij,i)) &
164     /(p(ij,1)-p(ij,jtt(ij)))
165     endif
166     endif
167     !
168     ! *** Find mixing ratio of precipitating downdraft ***
169     !
170     if(i.ne.inb(ij))then
171     if(i.eq.1)then
172     qstm=qs(ij,1)
173     else
174     qstm=qs(ij,i-1)
175     endif
176     if(mp(ij,i).gt.mp(ij,i+1))then
177     rat=mp(ij,i+1)/mp(ij,i)
178     qp(ij,i)=qp(ij,i+1)*rat+q(ij,i)*(1.0-rat)+100.*ginv* &
179     sigd*(ph(ij,i)-ph(ij,i+1))*(evap(ij,i)/mp(ij,i))
180     up(ij,i)=up(ij,i+1)*rat+u(ij,i)*(1.-rat)
181     vp(ij,i)=vp(ij,i+1)*rat+v(ij,i)*(1.-rat)
182     else
183     if(mp(ij,i+1).gt.0.0)then
184     qp(ij,i)=(gz(ij,i+1)-gz(ij,i) &
185     +qp(ij,i+1)*(lv(ij,i+1)+t(ij,i+1) &
186     *(cl-cpd))+cpd*(t(ij,i+1)-t(ij,i))) &
187     /(lv(ij,i)+t(ij,i)*(cl-cpd))
188     up(ij,i)=up(ij,i+1)
189     vp(ij,i)=vp(ij,i+1)
190     endif
191     endif
192     qp(ij,i)=min(qp(ij,i),qstm)
193     qp(ij,i)=max(qp(ij,i),0.0)
194     endif
195     endif
196     894 continue
197     899 continue
198     !
199     return
200     end

  ViewVC Help
Powered by ViewVC 1.1.21