module micromod  3
implicit none
contains


subroutine clean_up,5
use kvar; use kcnst
call pos_def(qrN,qvN,qvbank)
call pos_def(nrN,ccnN,ccnbank)
call pos_def(ncN,ccnN,ccnbank)
end subroutine clean_up


subroutine pos_def(q2,q1,qbank) 3,1
use kvar, only:mrl,mf,nu,k
real(mrl), dimension(:), intent(inout) :: q2,q1
real(mrl) :: qbank
	do k=1,nu
	if (q2(k).lt.0._mrl) then
		q1(k)=q1(k)+q2(k)
		q2(k)=0.
	end if
	if (q1(k).lt.0._mrl) then
		qbank=qbank+q1(k)*(mf(k+1)-mf(k))
		q1(k)=0.
	end if
	end do
end subroutine pos_def


subroutine microphysics(level) 3,7
use kvar; use kcnst
real(mrl) temx,presx,rhox,qsat
integer level
do k=1,nu
	temx=th(k)*exner(k)
	presx=psur*exner(k)**3.5
	rhox=(mf(k+1)-mf(k))/(zf(k+1)-zf(k))
	call sat_adj(temx,presx,qv(k),qc(k),ccn(k),nc(k),qsat)
	if (level.gt.1) then
		call driz_evap(temx,presx,rhox,qv(k),qr(k),nr(k),ccn(k),dt,qsat)
		call auto_conv(qr(k),nr(k),qc(k),nc(k),rhox,dt)
		call accretion(qr(k),qc(k),nc(k),dt)
		call driz_sed(qr(k),nr(k),qr_sed(k),nr_sed(k),droprad(k),rhox)
	end if
	th(k)=temx/exner(k)
	rh(k)=(qc(k)+qv(k)+qr(k))/qsat
 !       cloudfrac(k)=1./(1+exp(fracfac*(1.-rh(k))))
        cloudfrac(k)=min(1._mrl,max((rh(k)-fracfac)/(1-fracfac),0._mrl))
        if (fracfac.lt.0.) then
          cloudfrac(k)=0.
          if (qc(k).gt.1.e-12) cloudfrac(k)=1.
        endif

end do
do k=1,nu
  qtot(k)=qv(k)+qc(k)+qr(k)
  thetal(k)=(th(k)-elcp*(qc(k)+qr(k))/(exner(k+1)/2+exner(k)/2))*(1+.61*qv(k)-qc(k)-qr(k))
end do
end subroutine microphysics


subroutine driz_evap(t,p,r,qv,qr,nr,ccn,dt,qsat) 1,4
use precis; use kcnst; use kvar, only: elcp;
real(mrl) t,p,r,qv,qr,nr,ccn,dt,cc
real(mrl) esl,desdtl,qsat,dqrdt,c,dqr,dnr
call sat_water(t,p,esl,desdtl,qsat)
cc=(.5543903e11/t**2+2.042e7*t/(10*esl))
c=2.6388*(1.3333*pi*1.e3/r)**0.666667/cc
if (qv.lt.qsat.and.qr.lt.1.e-16) then
	qv=qv+qr
	qr=0.
	ccn=ccn+nr
	nr=0.
else
	dqrdt=c*(qv/qsat-1.)*(nr*nr*qr)**0.333333
	dqr=max(-qr,dqrdt*dt)
	dnr=nr*dqr/(qr+1.e-20)
	qv=qv-dqr
	t=t+elcp*dqr
	qr=qr+dqr
	nr=nr+dnr
	ccn=ccn-dnr
end if
end subroutine driz_evap 


subroutine accretion(qr,qc,nc,dt) 1,1
use precis
real(mrl), parameter:: mass_driz=6.545e-11
real(mrl) qr,qc,nc,dt
real(mrl) dqr
if (qc.gt.0.and.qr.gt.0.) then
	dqr=min(qc,dt*67.*(qc*qr)**1.15)
	nc=nc-nc*dqr/qc
	qr=qr+dqr
	qc=qc-dqr
end if
end subroutine accretion


subroutine driz_sed(qr,nr,qr_sed,nr_sed,rvr,rhox) 1,2
use precis; use kcnst
real(mrl), parameter:: mass_driz=6.545e-11
real(mrl) qr,nr,qr_sed,nr_sed,rhox,rvr
rvr=0.
if (nr.gt.0.) rvr=(3.*qr*rhox/(4.*pi*1.e3*nr*1.e6))**.33333
rvr=min(rvr,3.e-4_mrl)
qr_sed=-max(0._mrl,7000.*rvr-.1_mrl)
nr_sed=-max(0._mrl,12000.*rvr-.2_mrl)
end subroutine driz_sed 


subroutine auto_conv(qr,nr,qc,nc,rhox,dt) 1,1
use precis
real(mrl), parameter:: mass_driz=6.545e-11
real(mrl) qr,nr,qc,nc,rhox,dt
real(mrl) dqr,rvc,rvr
dqr=0.
if (qc.gt.0.and.nc.gt.0.) then
	dqr=min(qc,dt*1350*qc**2.47*nc**(-1.79))
	nr=nr+1.e-6*rhox*dqr/mass_driz
	nc=nc-nc*dqr/qc
	qr=qr+dqr
	qc=qc-dqr
end if
end subroutine auto_conv 


subroutine sat_adj(t,p,qv,ql,ccn,nc,qsat) 1,4
use precis; use kcnst; use kvar, only: elcp
real(mrl) t,p,qv,ql,dqv,qsat,desdtl,esl,ccn,nc,dnc
call sat_water(t,p,esl,desdtl,qsat)
	dqv=(qsat-qv)/(1+.611*desdtl*elcp/p)
	dqv=min(ql,dqv)
	if (dqv.lt.0.) then
		nc=nc+ccn
		ccn=0
	else
		dnc=-nc*dqv/(ql+1.e-20)
		nc=nc+dnc
		ccn=ccn-dnc
	end if
	qv=qv+dqv
	ql=ql-dqv
	t=t-elcp*dqv
if (ccn.lt.0.) ccn=0.
if (nc.lt.0.) nc=0.
if (qv.lt.0.) qv=0.
if (ql.lt.0.) ql=0.
	if (ql.lt.1.e-8) then
		ql=0.
		nc=0.
	end if	
end subroutine sat_adj


subroutine sat_water(t,p,esl,desdtl,qsatl) 3,1
use precis
real(mrl) t,p,esl,desdtl,qsatl
real(mrl) satfwa, satfwb
parameter ( satfwa = 1.0007 )
parameter ( satfwb = 3.46e-8 )  ! for p in Pa
real(mrl) satewa, satewb, satewc
parameter ( satewa = 611.21 )   ! es in Pa
parameter ( satewb = 17.502 )
parameter ( satewc = 32.18 )
!arps:
esl = (satfwa + satfwb * p) * satewa * exp( satewb*(t-273.15)/(t-satewc) )
desdtl = esl*satewb*(273.15-satewc)/(t-satewc)**2
!coamps:
!esl=611.0*10**(7.5*(t-273.16)/(t-35.86))
!desdtl = esl*2.302585*7.5*237.3/(t-35.86)**2
qsatl=.611*esl/p
return
end subroutine sat_water 


subroutine sat_ice(t,p,esi,desdti,qsati),1
use precis
real(mrl) t,p,esi,desdti,qsati
real(mrl) satfia, satfib
parameter ( satfia = 1.0003 )
parameter ( satfib = 4.18e-8 )  ! for p in Pa
real(mrl) sateia, sateib, sateic
parameter ( sateia = 611.15 )   ! es in Pa
parameter ( sateib = 22.452 )
parameter ( sateic = 0.6 )
esi = (satfia + satfib * p) * sateia * exp( sateib*(t-273.15)/(t-sateic) )
desdti = esi*sateib*(273.15-sateic)/(t-sateic)**2
qsati=.611*esi/p
return
end subroutine sat_ice

end module micromod