subroutine adprecpd( im, ix, km, dt, del, sl, ps, rhc, q_in, & 2,4
           cwm_in, t_in, q_out, cwm_out, t_out, rn_out, adq_in, adcwm_in, &
           adt_in, adq_out, adcwm_out, adt_out, adrn_out, &
           adjoint,mype)
!C***************************************************************
!C***************************************************************
!C** This routine was generated by the                         **
!C** Tangent linear and Adjoint Model Compiler,  TAMC 5.3.0    **
!C***************************************************************
!C***************************************************************
!C==============================================
!C all entries are defined explicitly
!C==============================================
 use type_kinds, only: fp_kind
      implicit none
      include 'constant.h'
!C==============================================
!C define parameters                            
!C==============================================
      logical adjoint0
      parameter (adjoint0=.false.)
      real(fp_kind) c1
      parameter ( c1 = 300. )
      real(fp_kind) c2
      parameter ( c2 = 0.5 )
      real(fp_kind) climit
      parameter ( climit = 1.e-20 )
      real(fp_kind) cmr
      parameter ( cmr = 1./0.0003 )
      real(fp_kind) cws
      parameter ( cws = 0.025 )
      real(fp_kind) d00
      parameter ( d00 = 0. )
      real(fp_kind) eliv
      parameter ( eliv = hvap+hfus )
      real(fp_kind) elwv
      parameter ( elwv = hvap )
      real(fp_kind) eliw
      parameter ( eliw = eliv-elwv )
      real(fp_kind) eps
      parameter ( eps = rd/rv )
      real(fp_kind) epsm1
      parameter ( epsm1 = eps-1. )
      real(fp_kind) epsq
      parameter ( epsq = 2.e-12 )
      real(fp_kind) g
      parameter ( g = grav )
      real(fp_kind) h1
      parameter ( h1 = 1. )
      real(fp_kind) h1000
      parameter ( h1000 = 1000. )
      real(fp_kind) ke
      parameter ( ke = 0.00002 )
      real(fp_kind) rcp
      parameter ( rcp = h1/cp )
      real(fp_kind) row
      parameter ( row = 1000. )
      real(fp_kind) rrow
      parameter ( rrow = h1/row )
      real(fp_kind) us
      parameter ( us = h1 )

!C==============================================
!C define arguments
!C==============================================
      logical adjoint
      integer mype
      integer ix
      integer km
      real(fp_kind) adcwm_in(km,ix)
      real(fp_kind) adcwm_out(km,ix)
      real(fp_kind) adq_in(km,ix)
      real(fp_kind) adq_out(km,ix)
      integer im
      real(fp_kind) adrn_out(im)
      real(fp_kind) adt_in(km,ix)
      real(fp_kind) adt_out(km,ix)
      real(fp_kind) cwm_in(km,ix)
      real(fp_kind) cwm_out(km,ix)
      real(fp_kind) del(km)
      real(fp_kind) dt
      real(fp_kind) ps(im)
      real(fp_kind) q_in(km,ix)
      real(fp_kind) q_out(km,ix)
      real(fp_kind) rhc(km)
      real(fp_kind) rn_out(im)
      real(fp_kind) sl(km)
      real(fp_kind) t_in(km,ix)
      real(fp_kind) t_out(km,ix)

!C==============================================
!C define local variables
!C==============================================
      real(fp_kind) aa2
      real(fp_kind) adamaxcmr
      real(fp_kind) adamaxcms
      real(fp_kind) adamaxps
      real(fp_kind) adamaxrq
      real(fp_kind) adccr
      real(fp_kind) adcwmin
      real(fp_kind) adcwmk
      real(fp_kind) adcwmks
      real(fp_kind) aderk
      real(fp_kind) aderr
      real(fp_kind) aders
      real(fp_kind) ades
      real(fp_kind) adexpf
      real(fp_kind) adfi
      real(fp_kind) adppr0
      real(fp_kind) adppr1
      real(fp_kind) adppr2
      real(fp_kind) adpps0
      real(fp_kind) adpps1
      real(fp_kind) adpps2
      real(fp_kind) adpraut
      real(fp_kind) adpraut0
      real(fp_kind) adprecrk0
      real(fp_kind) adprecrk1
      real(fp_kind) adprecrl1k
      real(fp_kind) adprecrl_0
      real(fp_kind) adprecrl_1
      real(fp_kind) adprecrl_2
      real(fp_kind) adprecrl_3
      real(fp_kind) adprecsk0
      real(fp_kind) adprecsk1
      real(fp_kind) adprecsl1k
      real(fp_kind) adprecsl_0
      real(fp_kind) adprecsl_1
      real(fp_kind) adprecsl_2
      real(fp_kind) adprecsl_3
      real(fp_kind) adpsaci
      real(fp_kind) adpsaut
      real(fp_kind) adpsm1
      real(fp_kind) adpsm2
      real(fp_kind) adpsm3
      real(fp_kind) adq2
      real(fp_kind) adqin
      real(fp_kind) adqk
      real(fp_kind) adqs
      real(fp_kind) adqs0
      real(fp_kind) adrprs
      real(fp_kind) adrq
      real(fp_kind) adrqkll
      real(fp_kind) adt2
      real(fp_kind) adtem1
      real(fp_kind) adtem2
      real(fp_kind) adtem3
      real(fp_kind) adtem4
      real(fp_kind) adtem5
      real(fp_kind) adterm
      real(fp_kind) adterm1
      real(fp_kind) adterm2r
      real(fp_kind) adterm2s
      real(fp_kind) adterm3r
      real(fp_kind) adterm3s
      real(fp_kind) adterm4
      real(fp_kind) adterm5
      real(fp_kind) adterm6
      real(fp_kind) adtin
      real(fp_kind) adtmt0
      real(fp_kind) adtmt0k
      real(fp_kind) adww0
      real(fp_kind) adww1
      real(fp_kind) adww2
      real(fp_kind) adwwn
      real(fp_kind) amaxcmr
      real(fp_kind) amaxcms
      real(fp_kind) amaxps
      real(fp_kind) amaxrq
      real(fp_kind) c00
      real(fp_kind) ccr
      logical comput,comput0(ix)
      real(fp_kind) conde
      real(fp_kind) condt
      real(fp_kind) const
      real(fp_kind) cr
      real(fp_kind) crs1
      real(fp_kind) crs2
      real(fp_kind) csm1
      real(fp_kind) cwmin
      real(fp_kind) cwmk
      real(fp_kind) cwmks
      real(fp_kind) dtcp
      real(fp_kind) dum1
      real(fp_kind) dum2
      real(fp_kind) erk
      real(fp_kind) err
      real(fp_kind) ers
      real(fp_kind) es
      real(fp_kind) expf
      real(fp_kind) fi
      integer i
      integer iwl
      integer iwl1,iwl1k(km+1,ix)
      integer k
      integer k2
      real(fp_kind) ppr0
      real(fp_kind) ppr1
      real(fp_kind) ppr2
      real(fp_kind) pps0
      real(fp_kind) pps1
      real(fp_kind) pps2
      real(fp_kind) praut
      real(fp_kind) praut0
      real(fp_kind) precrk0
      real(fp_kind) precrk1
      real(fp_kind) precrl1k,precrl1ki(km+1,ix)
      real(fp_kind) precrl_0
      real(fp_kind) precrl_1
      real(fp_kind) precrl_2
      real(fp_kind) precrl_3
      real(fp_kind) precsk0
      real(fp_kind) precsk1
      real(fp_kind) precsl1k,precsl1ki(km+1,ix)
      real(fp_kind) precsl_0
      real(fp_kind) precsl_1
      real(fp_kind) precsl_2
      real(fp_kind) precsl_3
      real(fp_kind) pres
      real(fp_kind) psaci
      real(fp_kind) psaut
      real(fp_kind) psm1
      real(fp_kind) psm2
      real(fp_kind) psm3
      real(fp_kind) q2
      real(fp_kind) qin
      real(fp_kind) qk
      real(fp_kind) qs
      real(fp_kind) qs0
      real(fp_kind) rconde
      real(fp_kind) rdt
      real(fp_kind) rhci(km)
      real(fp_kind) rke
      real(fp_kind) rprs
      real(fp_kind) rq
      real(fp_kind) rqkll
      real(fp_kind) t2
      real(fp_kind) tem
      real(fp_kind) tem1
      real(fp_kind) tem2
      real(fp_kind) tem3
      real(fp_kind) tem4
      real(fp_kind) tem5
      real(fp_kind) term
      real(fp_kind) term1
      real(fp_kind) term2r
      real(fp_kind) term2s
      real(fp_kind) term3r
      real(fp_kind) term3s
      real(fp_kind) term4
      real(fp_kind) term5
      real(fp_kind) term6
      real(fp_kind) tin
      real(fp_kind) tmt0
      real(fp_kind) tmt0k
      real(fp_kind) wmin
      real(fp_kind) ww0
      real(fp_kind) ww1
      real(fp_kind) ww2
      real(fp_kind) wwn
      real(fp_kind) zaodt

!C----------------------------------------------
!C RESET LOCAL ADJOINT VARIABLES
!C----------------------------------------------
      adamaxcmr = 0.
      adamaxcms = 0.
      adamaxps = 0.
      adamaxrq = 0.
      adccr = 0.
      adcwmin = 0.
      adcwmk = 0.
      adcwmks = 0.
      aderk = 0.
      aderr = 0.
      aders = 0.
      ades = 0.
      adexpf = 0.
      adfi = 0.
      adppr0 = 0.
      adppr1 = 0.
      adppr2 = 0.
      adpps0 = 0.
      adpps1 = 0.
      adpps2 = 0.
      adpraut = 0.
      adpraut0 = 0.
      adprecrk0 = 0.
      adprecrk1 = 0.
      adprecrl1k = 0.
      adprecrl_0 = 0.
      adprecrl_1 = 0.
      adprecrl_2 = 0.
      adprecrl_3 = 0.
      adprecsk0 = 0.
      adprecsk1 = 0.
      adprecsl1k = 0.
      adprecsl_0 = 0.
      adprecsl_1 = 0.
      adprecsl_2 = 0.
      adprecsl_3 = 0.
      adpsaci = 0.
      adpsaut = 0.
      adpsm1 = 0.
      adpsm2 = 0.
      adpsm3 = 0.
      adq2 = 0.
      adqin = 0.
      adqk = 0.
      adqs = 0.
      adqs0 = 0.
      adrprs = 0.
      adrq = 0.
      adrqkll = 0.
      adt2 = 0.
      adtem1 = 0.
      adtem2 = 0.
      adtem3 = 0.
      adtem4 = 0.
      adtem5 = 0.
      adterm = 0.
      adterm1 = 0.
      adterm2r = 0.
      adterm2s = 0.
      adterm3r = 0.
      adterm3s = 0.
      adterm4 = 0.
      adterm5 = 0.
      adterm6 = 0.
      adtin = 0.
      adtmt0 = 0.
      adtmt0k = 0.
      adww0 = 0.
      adww1 = 0.
      adww2 = 0.
      adwwn = 0.

!C----------------------------------------------
!C ROUTINE BODY
!C----------------------------------------------
!C----------------------------------------------
!C FUNCTION AND TAPE COMPUTATIONS
!C----------------------------------------------
      rdt = h1/dt
      zaodt = 800.*rdt
      csm1 = 5.e-8*zaodt
      crs1 = 5.e-6*zaodt
      crs2 = 6.666e-10*zaodt
      cr = 0.0005*zaodt
      aa2 = 0.00125*zaodt
      rke = ke*sqrt(rdt)
      dtcp = dt*rcp
      c00 = 0.0001*dt
      do k = 1, km
        rhci(k) = 1./rhc(k)
      end do
      do i = 1, im

!!new
        comput0(i) =  .false.
        do k = 1, km
          tem = 0.00001*sl(k)*ps(i)*0.01
          if (cwm_in(k,i) .gt. tem) then
            comput0(i) =  .true.
          endif
        end do
!!new
         
        iwl1 = 0
        precrl1k = 0.
        precsl1k = 0.
        rn_out(i) = 0.
        const = ps(i)*(h1000*dt/g)

!!new
        if (comput0(i)) then
!!new
        do k = km, 1, -1
          wmin = 0.00001*sl(k)*ps(i)*0.01
          precrl_0 = precrl1k
          precsl_0 = precsl1k

!         Save forward model values for use in adjoint model
          precrl1ki(k,i) = precrl_0
          precsl1ki(k,i) = precsl_0
          iwl1k(k,i)     = iwl1

          iwl = 0
          tin = t_in(k,i)
          qin = q_in(k,i)
          cwmin = cwm_in(k,i)
          pres = ps(i)*sl(k)
          if (precrl_0 .gt. 0.) then
            precrk0 = precrl_0
          else
            precrk0 = 0.
          endif
          if (precsl_0 .gt. 0.) then
            precsk0 = precsl_0
          else
            precsk0 = 0.
          endif
          if (cwmin .gt. climit) then
            wwn = cwmin
          else
            wwn = climit
          endif
          if (wwn .gt. climit .or. precrk0+precsk0 .gt. d00) then
            comput =  .true. 
          else
            comput =  .false. 
          endif
          if (comput) then
            conde = const*del(k)
            condt = conde*rdt
            rconde = h1/conde
            if (qin .gt. epsq) then
              qk = qin
            else
              qk = epsq
            endif
            tmt0 = tin-273.16
            call adfpvsx(tin,es,dum1,dum2,adjoint0)
            qs0 = eps*es/(pres+epsm1*es)
            if (qs0 .gt. epsq) then
              qs = qs0
            else
              qs = epsq
            endif
            if (tmt0 .lt. (-15.)) then
              fi = qk-rhc(k)*qs
              if (fi .gt. d00 .or. wwn .gt. climit) then
                iwl = 1
              else
                iwl = 0
              endif
            else if (tmt0 .ge. 0.) then
              iwl = 0
            else
              iwl = 0
              if (iwl1 .eq. 1 .and. wwn .gt. climit) then
                iwl = 1
              endif
            endif
            if (qs .le. 1.e-10) then
              rq = d00
            else
              rq = qk/qs
            endif
            if (rq .lt. rhc(k)) then
              ccr = d00
            else if (rq .ge. us) then
              ccr = us
            else
              if (rq .lt. us) then
                rqkll = rq
              else
                rqkll = us
              endif
              ccr = h1-sqrt((us-rqkll)/(us-rhc(k)))
            endif
            if (ccr .gt. 0.) then
              ww0 = cwmin
              if (ww0 .gt. 0.) then
                cwmk = ww0
              else
                cwmk = 0.
              endif
              if (iwl .eq. 1) then
                term1 = cwmk-wmin
                if (term1 .gt. 0.) then
                  amaxcms = term1
                else
                  amaxcms = 0.
                endif
                expf = dt*exp(0.025*tmt0)
                term2s = 0.0004*expf*amaxcms
                if (term2s .lt. cwmk) then
                  psaut = term2s
                else
                  psaut = cwmk
                endif
                ww1 = ww0-psaut
                if (ww1 .gt. 0.) then
                  cwmks = ww1
                else
                  cwmks = 0.
                endif
                term3s = aa2*expf*precsl_0*cwmks
                if (term3s .lt. cwmks) then
                  psaci = term3s
                else
                  psaci = cwmks
                endif
                ww2 = ww1-psaci
                precsl_1 = precsl_0+(ww0-ww2)*condt
                precrl_1 = precrl_0
              else
                amaxcmr = cwmk
                tem1 = precsl_0+precrl_0
                term2r = 268.-tin
                if (term2r .gt. 0.) then
                  term3r = term2r
                else
                  term3r = 0.
                endif
                if (term3r .lt. 20.) then
                  term4 = term3r
                else
                  term4 = 20.
                endif
                tem2 = term4
                tem3 = (1.+c1*sqrt(tem1*rdt))*(1+c2*sqrt(tem2))
                if (ccr .gt. 0.01) then
                  term5 = ccr
                else
                  term5 = 0.01
                endif
                tem4 = amaxcmr*cmr*tem3/term5
                term6 = tem4*tem4
                if (term6 .lt. 50.) then
                  tem5 = term6
                else
                  tem5 = 50.
                endif
                praut = c00*tem3*amaxcmr*(1.-exp(-tem5))
                if (praut .lt. cwmk) then
                  praut0 = praut
                else
                  praut0 = cwmk
                endif
                ww2 = ww0-praut0
                precrl_1 = precrl_0+(ww0-ww2)*condt
                precsl_1 = precsl_0
              endif
            else
              ww2 = cwmin
              precrl_1 = precrl_0
              precsl_1 = precsl_0
            endif
            if (tmt0 .gt. (-30.)) then
              tmt0k = tmt0
            else
              tmt0k = -30.
            endif
            if (precrl_1 .gt. 0.) then
              precrk1 = precrl_1
            else
              precrk1 = 0.
            endif
            if (precsl_1 .gt. 0.) then
              precsk1 = precsl_1
            else
              precsk1 = 0.
            endif
            term = rhc(k)-rq
            if (term .gt. 0.) then
              amaxrq = term*conde
            else
              amaxrq = 0.
            endif
            ppr0 = rke*amaxrq*sqrt(precrk1)
            if (tmt0 .ge. 0.) then
              pps0 = 0.
            else
              pps0 = (crs1+crs2*tmt0k)*amaxrq*precsk1*rhci(k)
            endif
            erk = precrk1+precsk1
            if (rq .ge. 1.e-10) then
              erk = amaxrq*qk*rdt/rq
            endif
            if (ppr0+pps0 .gt. abs(erk)) then
              rprs = erk/(precrk1+precsk1)
              ppr1 = precrk1*rprs
              pps1 = precsk1*rprs
            else
              ppr1 = ppr0
              pps1 = pps0
            endif
            if (ppr1 .lt. precrk1) then
              ppr2 = ppr1
            else
              ppr2 = precrk1
            endif
            if (pps1 .lt. precsk1) then
              pps2 = pps1
            else
              pps2 = precsk1
            endif
            err = ppr2*rconde
            ers = pps2*rconde
            precrl_2 = precrl_1-ppr2
            precsl_2 = precsl_1-pps2
            if (tmt0 .gt. 0.) then
              if (precsl_2 .gt. 0.) then
                amaxps = precsl_2
              else
                amaxps = 0.
              endif
              psm1 = csm1*tmt0*tmt0*amaxps
              if (ww2 .gt. 0.) then
                psm2 = cws*cr*ww2*amaxps
              else
                psm2 = 0.
              endif
              ppr0 = (psm1+psm2)*conde
              if (ppr0 .gt. amaxps) then
                ppr1 = amaxps
                psm3 = amaxps*rconde
              else
                ppr1 = ppr0
                psm3 = psm1
              endif
              precrl_3 = precrl_2+ppr1
              precsl_3 = precsl_2-ppr1
            else
              psm3 = d00
              precrl_3 = precrl_2
              precsl_3 = precsl_2
            endif
            t2 = tin-dtcp*(elwv*err+eliv*ers+eliw*psm3)
            q2 = qin+dt*(err+ers)
          else
            t2 = tin
            q2 = qin
            precrl_3 = precrl_0
            precsl_3 = precsl_0
            ww2 = cwmin
          endif
          iwl1 = iwl
          if (precrl_3 .gt. 0.) then
            precrl1k = precrl_3
          else
            precrl1k = 0.
          endif
          if (precsl_3 .gt. 0.) then
            precsl1k = precsl_3
          else
            precsl1k = 0.
          endif
          if (ww2 .lt. 0.) then
            q_out(k,i) = q2+ww2
            t_out(k,i) = t2-elwv*rcp*ww2
            cwm_out(k,i) = 0.
          else
            q_out(k,i) = q2
            t_out(k,i) = t2
            cwm_out(k,i) = ww2
          endif
        end do
        rn_out(i) = (precrl1k+precsl1k)*rrow

!!new - no compute case
      else
         do k = km, 1, -1
            if (cwm_in(k,i) .lt. 0.) then
              q_out(k,i)   = q_in(k,i)+cwm_in(k,i)
              t_out(k,i)   = t_in(k,i)-elwv*rcp*cwm_in(k,i)
              cwm_out(k,i) = 0.
            else
              q_out(k,i)   = q_in(k,i)
              t_out(k,i)   = t_in(k,i)
              cwm_out(k,i) = cwm_in(k,i)
            endif
          end do
          rn_out(i) = 0.
        endif

!     End of loop over profiles         
      end do

      if (.not.adjoint) return

!C----------------------------------------------
!C ADJOINT COMPUTATIONS
!C----------------------------------------------
      do i = im, 1, -1

       if (comput0(i)) then

        const = ps(i)*(h1000*dt/g)
        adprecrl1k = adprecrl1k+adrn_out(i)*rrow
        adprecsl1k = adprecsl1k+adrn_out(i)*rrow
        adrn_out(i) = 0.
        do k = 1, km

!         Load saved forward model values
          precrl1k = precrl1ki(k,i)
          precsl1k = precsl1ki(k,i)
          iwl1     = iwl1k(k,i)

!         Recompute forward model values
          wmin = 0.00001*sl(k)*ps(i)*0.01
          precrl_0 = precrl1k
          precsl_0 = precsl1k
          tin = t_in(k,i)
          qin = q_in(k,i)
          cwmin = cwm_in(k,i)
          pres = ps(i)*sl(k)
          if (precrl_0 .gt. 0.) then
            precrk0 = precrl_0
          else
            precrk0 = 0.
          endif
          if (precsl_0 .gt. 0.) then
            precsk0 = precsl_0
          else
            precsk0 = 0.
          endif
          if (cwmin .gt. climit) then
            wwn = cwmin
          else
            wwn = climit
          endif
          if (wwn .gt. climit .or. precrk0+precsk0 .gt. d00) then
            comput =  .true. 
          else
            comput =  .false. 
          endif
          if (comput) then
            conde = const*del(k)
            condt = conde*rdt
            if (qin .gt. epsq) then
              qk = qin
            else
              qk = epsq
            endif
            tmt0 = tin-273.16
            call adfpvsx(tin,es,dum1,dum2,adjoint0)
            qs0 = eps*es/(pres+epsm1*es)
            if (qs0 .gt. epsq) then
              qs = qs0
            else
              qs = epsq
            endif
            if (tmt0 .lt. (-15.)) then
              fi = qk-rhc(k)*qs
              if (fi .gt. d00 .or. wwn .gt. climit) then
                iwl = 1
              else
                iwl = 0
              endif
            else if (tmt0 .ge. 0.) then
              iwl = 0
            else
              iwl = 0
              if (iwl1 .eq. 1 .and. wwn .gt. climit) then
                iwl = 1
              endif
            endif
            if (qs .le. 1.e-10) then
              rq = d00
            else
              rq = qk/qs
            endif
            if (rq .lt. rhc(k)) then
              ccr = d00
            else if (rq .ge. us) then
              ccr = us
            else
              if (rq .lt. us) then
                rqkll = rq
              else
                rqkll = us
              endif
              ccr = h1-sqrt((us-rqkll)/(us-rhc(k)))
            endif
            if (ccr .gt. 0.) then
              ww0 = cwmin
              if (ww0 .gt. 0.) then
                cwmk = ww0
              else
                cwmk = 0.
              endif
              if (iwl .eq. 1) then
                term1 = cwmk-wmin
                if (term1 .gt. 0.) then
                  amaxcms = term1
                else
                  amaxcms = 0.
                endif
                expf = dt*exp(0.025*tmt0)
                term2s = 0.0004*expf*amaxcms
                if (term2s .lt. cwmk) then
                  psaut = term2s
                else
                  psaut = cwmk
                endif
                ww1 = ww0-psaut
                if (ww1 .gt. 0.) then
                  cwmks = ww1
                else
                  cwmks = 0.
                endif
                term3s = aa2*expf*precsl_0*cwmks
                if (term3s .lt. cwmks) then
                  psaci = term3s
                else
                  psaci = cwmks
                endif
                ww2 = ww1-psaci
                precsl_1 = precsl_0+(ww0-ww2)*condt
                precrl_1 = precrl_0
              else
                amaxcmr = cwmk
                tem1 = precsl_0+precrl_0
                term2r = 268.-tin
                if (term2r .gt. 0.) then
                  term3r = term2r
                else
                  term3r = 0.
                endif
                if (term3r .lt. 20.) then
                  term4 = term3r
                else
                  term4 = 20.
                endif
                tem2 = term4
                tem3 = (1.+c1*sqrt(tem1*rdt))*(1+c2*sqrt(tem2))
                if (ccr .gt. 0.01) then
                  term5 = ccr
                else
                  term5 = 0.01
                endif
                tem4 = amaxcmr*cmr*tem3/term5
                term6 = tem4*tem4
                if (term6 .lt. 50.) then
                  tem5 = term6
                else
                  tem5 = 50.
                endif
                praut = c00*tem3*amaxcmr*(1.-exp(-tem5))
                if (praut .lt. cwmk) then
                  praut0 = praut
                else
                  praut0 = cwmk
                endif
                ww2 = ww0-praut0
                precrl_1 = precrl_0+(ww0-ww2)*condt
                precsl_1 = precsl_0
              endif
            else
              ww2 = cwmin
              precrl_1 = precrl_0
              precsl_1 = precsl_0
            endif
            if (tmt0 .gt. (-30.)) then
              tmt0k = tmt0
            else
              tmt0k = -30.
            endif
            if (precrl_1 .gt. 0.) then
              precrk1 = precrl_1
            else
              precrk1 = 0.
            endif
            if (precsl_1 .gt. 0.) then
              precsk1 = precsl_1
            else
              precsk1 = 0.
            endif
            term = rhc(k)-rq
            if (term .gt. 0.) then
              amaxrq = term*conde
            else
              amaxrq = 0.
            endif
!!            ppr0 = rke*amaxrq*sqrt(precrk1)
            if (precrk1>0) then
               ppr0 = rke*amaxrq*sqrt(precrk1)
            else
               ppr0 = 0.0
            endif

            if (tmt0 .ge. 0.) then
              pps0 = 0.
            else
              pps0 = (crs1+crs2*tmt0k)*amaxrq*precsk1*rhci(k)
            endif
            erk = precrk1+precsk1
            if (rq .ge. 1.e-10) then
              erk = amaxrq*qk*rdt/rq
            endif
            if (ppr0+pps0 .gt. abs(erk)) then
              rprs = erk/(precrk1+precsk1)
              ppr1 = precrk1*rprs
              pps1 = precsk1*rprs
            else
              ppr1 = ppr0
              pps1 = pps0
            endif
            if (ppr1 .lt. precrk1) then
              ppr2 = ppr1
            else
              ppr2 = precrk1
            endif
            if (pps1 .lt. precsk1) then
              pps2 = pps1
            else
              pps2 = precsk1
            endif
            precrl_2 = precrl_1-ppr2
            precsl_2 = precsl_1-pps2
            if (tmt0 .gt. 0.) then
              if (precsl_2 .gt. 0.) then
                amaxps = precsl_2
              else
                amaxps = 0.
              endif
              psm1 = csm1*tmt0*tmt0*amaxps
              if (ww2 .gt. 0.) then
                psm2 = cws*cr*ww2*amaxps
              else
                psm2 = 0.
              endif
              ppr0 = (psm1+psm2)*conde
              if (ppr0 .gt. amaxps) then
                ppr1 = amaxps
              else
                ppr1 = ppr0
              endif
              precrl_3 = precrl_2+ppr1
              precsl_3 = precsl_2-ppr1
            else
              precrl_3 = precrl_2
              precsl_3 = precsl_2
            endif
          else
            precrl_3 = precrl_0
            precsl_3 = precsl_0
            ww2 = cwmin
          endif

!         Adjoint model
          if (ww2 .lt. 0.) then
            adcwm_out(k,i) = 0.
            adt2 = adt2+adt_out(k,i)
            adww2 = adww2-adt_out(k,i)*elwv*rcp
            adt_out(k,i) = 0.
            adq2 = adq2+adq_out(k,i)
            adww2 = adww2+adq_out(k,i)
            adq_out(k,i) = 0.
          else
            adww2 = adww2+adcwm_out(k,i)
            adcwm_out(k,i) = 0.
            adt2 = adt2+adt_out(k,i)
            adt_out(k,i) = 0.
            adq2 = adq2+adq_out(k,i)
            adq_out(k,i) = 0.
          endif
          if (precsl_3 .gt. 0.) then
            adprecsl_3 = adprecsl_3+adprecsl1k
            adprecsl1k = 0.
          else
            adprecsl1k = 0.
          endif
          if (precrl_3 .gt. 0.) then
            adprecrl_3 = adprecrl_3+adprecrl1k
            adprecrl1k = 0.
          else
            adprecrl1k = 0.
          endif
          if (comput) then
            rconde = h1/conde
            aderr = aderr+adq2*dt
            aders = aders+adq2*dt
            adqin = adqin+adq2
            adq2 = 0.
            aderr = aderr-adt2*dtcp*elwv
            aders = aders-adt2*dtcp*eliv
            adpsm3 = adpsm3-adt2*dtcp*eliw
            adtin = adtin+adt2
            adt2 = 0.
            if (tmt0 .gt. 0.) then
              adppr1 = adppr1-adprecsl_3
              adprecsl_2 = adprecsl_2+adprecsl_3
              adprecsl_3 = 0.
              adppr1 = adppr1+adprecrl_3
              adprecrl_2 = adprecrl_2+adprecrl_3
              adprecrl_3 = 0.
              if (ppr0 .gt. amaxps) then
                adamaxps = adamaxps+adpsm3*rconde
                adpsm3 = 0.
                adamaxps = adamaxps+adppr1
                adppr1 = 0.
              else
                adpsm1 = adpsm1+adpsm3
                adpsm3 = 0.
                adppr0 = adppr0+adppr1
                adppr1 = 0.
              endif
              adpsm1 = adpsm1+adppr0*conde
              adpsm2 = adpsm2+adppr0*conde
              adppr0 = 0.
              if (ww2 .gt. 0.) then
                adamaxps = adamaxps+adpsm2*cws*cr*ww2
                adww2 = adww2+adpsm2*cws*cr*amaxps
                adpsm2 = 0.
              else
                adpsm2 = 0.
              endif
              adamaxps = adamaxps+adpsm1*csm1*tmt0*tmt0
              adtmt0 = adtmt0+2*adpsm1*csm1*tmt0*amaxps
              adpsm1 = 0.
              if (precsl_2 .gt. 0.) then
                adprecsl_2 = adprecsl_2+adamaxps
                adamaxps = 0.
              else
                adamaxps = 0.
              endif
            else
              adprecsl_2 = adprecsl_2+adprecsl_3
              adprecsl_3 = 0.
              adprecrl_2 = adprecrl_2+adprecrl_3
              adprecrl_3 = 0.
              adpsm3 = 0.
            endif
            adpps2 = adpps2-adprecsl_2
            adprecsl_1 = adprecsl_1+adprecsl_2
            adprecsl_2 = 0.
            adppr2 = adppr2-adprecrl_2
            adprecrl_1 = adprecrl_1+adprecrl_2
            adprecrl_2 = 0.
            adpps2 = adpps2+aders*rconde
            aders = 0.
            adppr2 = adppr2+aderr*rconde
            aderr = 0.
            if (pps1 .lt. precsk1) then
              adpps1 = adpps1+adpps2
              adpps2 = 0.
            else
              adprecsk1 = adprecsk1+adpps2
              adpps2 = 0.
            endif
            ppr0 = rke*amaxrq*sqrt(precrk1)
            if (ppr0+pps0 .gt. abs(erk)) then
              ppr1 = precrk1*rprs
            else
              ppr1 = ppr0
            endif
            if (ppr1 .lt. precrk1) then
              adppr1 = adppr1+adppr2
              adppr2 = 0.
            else
              adprecrk1 = adprecrk1+adppr2
              adppr2 = 0.
            endif
            ppr0 = rke*amaxrq*sqrt(precrk1)
            if (ppr0+pps0 .gt. abs(erk)) then
              adprecsk1 = adprecsk1+adpps1*rprs
              adrprs = adrprs+adpps1*precsk1
              adpps1 = 0.
              adprecrk1 = adprecrk1+adppr1*rprs
              adrprs = adrprs+adppr1*precrk1
              adppr1 = 0.
              aderk = aderk+adrprs/(precrk1+precsk1)
              adprecrk1 = adprecrk1-adrprs*(erk/((precrk1+precsk1)* &
                   (precrk1+precsk1)))
              adprecsk1 = adprecsk1-adrprs*(erk/((precrk1+precsk1)* &
                   (precrk1+precsk1)))
              adrprs = 0.
            else
              adpps0 = adpps0+adpps1
              adpps1 = 0.
              adppr0 = adppr0+adppr1
              adppr1 = 0.
            endif
            if (rq .ge. 1.e-10) then
              adamaxrq = adamaxrq+aderk*(qk*rdt/rq)
              adqk = adqk+aderk*(amaxrq*rdt/rq)
              adrq = adrq-aderk*(amaxrq*qk*rdt/(rq*rq))
              aderk = 0.
            endif
            adprecrk1 = adprecrk1+aderk
            adprecsk1 = adprecsk1+aderk
            aderk = 0.
            if (tmt0 .ge. 0.) then
              adpps0 = 0.
            else
              adamaxrq = adamaxrq+adpps0*(crs1+crs2*tmt0k)*precsk1* &
                   rhci(k)
              adprecsk1 = adprecsk1+adpps0*(crs1+crs2*tmt0k)*amaxrq* &
                   rhci(k)
              adtmt0k = adtmt0k+adpps0*crs2*amaxrq*precsk1*rhci(k)
              adpps0 = 0.
            endif
!!            adamaxrq = adamaxrq+adppr0*rke*sqrt(precrk1)
!!            adprecrk1 = adprecrk1+adppr0*rke*amaxrq*1./(2.* &
!!                 sqrt(precrk1))
!!            adppr0 = 0.
            if (precrk1>0) then
               adamaxrq = adamaxrq+adppr0*rke*sqrt(precrk1)
               adprecrk1 = adprecrk1+adppr0*rke*amaxrq*1./(2.* &
                    sqrt(precrk1))
               adppr0 = 0.
            else
               adppr0 = 0.
            endif

            if (term .gt. 0.) then
              adterm = adterm+adamaxrq*conde
              adamaxrq = 0.
            else
              adamaxrq = 0.
            endif
            adrq = adrq-adterm
            adterm = 0.
            if (precsl_1 .gt. 0.) then
              adprecsl_1 = adprecsl_1+adprecsk1
              adprecsk1 = 0.
            else
              adprecsk1 = 0.
            endif
            if (precrl_1 .gt. 0.) then
              adprecrl_1 = adprecrl_1+adprecrk1
              adprecrk1 = 0.
            else
              adprecrk1 = 0.
            endif
            if (tmt0 .gt. (-30.)) then
              adtmt0 = adtmt0+adtmt0k
              adtmt0k = 0.
            else
              adtmt0k = 0.
            endif
            if (ccr .gt. 0.) then
              if (iwl .eq. 1) then
                adprecrl_0 = adprecrl_0+adprecrl_1
                adprecrl_1 = 0.
                adprecsl_0 = adprecsl_0+adprecsl_1
                adww0 = adww0+adprecsl_1*condt
                adww2 = adww2-adprecsl_1*condt
                adprecsl_1 = 0.
                adpsaci = adpsaci-adww2
                adww1 = adww1+adww2
                adww2 = 0.
                if (term3s .lt. cwmks) then
                  adterm3s = adterm3s+adpsaci
                  adpsaci = 0.
                else
                  adcwmks = adcwmks+adpsaci
                  adpsaci = 0.
                endif
                adcwmks = adcwmks+adterm3s*aa2*expf*precsl_0
                adexpf = adexpf+adterm3s*aa2*precsl_0*cwmks
                adprecsl_0 = adprecsl_0+adterm3s*aa2*expf*cwmks
                adterm3s = 0.
                if (ww1 .gt. 0.) then
                  adww1 = adww1+adcwmks
                  adcwmks = 0.
                else
                  adcwmks = 0.
                endif
                adpsaut = adpsaut-adww1
                adww0 = adww0+adww1
                adww1 = 0.
                if (term2s .lt. cwmk) then
                  adterm2s = adterm2s+adpsaut
                  adpsaut = 0.
                else
                  adcwmk = adcwmk+adpsaut
                  adpsaut = 0.
                endif
                adamaxcms = adamaxcms+0.0004*adterm2s*expf
                adexpf = adexpf+0.0004*adterm2s*amaxcms
                adterm2s = 0.
                adtmt0 = adtmt0+0.025*adexpf*dt*exp(0.025*tmt0)
                adexpf = 0.
                if (term1 .gt. 0.) then
                  adterm1 = adterm1+adamaxcms
                  adamaxcms = 0.
                else
                  adamaxcms = 0.
                endif
                adcwmk = adcwmk+adterm1
                adterm1 = 0.
              else
                adprecsl_0 = adprecsl_0+adprecsl_1
                adprecsl_1 = 0.
                adprecrl_0 = adprecrl_0+adprecrl_1
                adww0 = adww0+adprecrl_1*condt
                adww2 = adww2-adprecrl_1*condt
                adprecrl_1 = 0.
                adpraut0 = adpraut0-adww2
                adww0 = adww0+adww2
                adww2 = 0.
                if (praut .lt. cwmk) then
                  adpraut = adpraut+adpraut0
                  adpraut0 = 0.
                else
                  adcwmk = adcwmk+adpraut0
                  adpraut0 = 0.
                endif
                adamaxcmr = adamaxcmr+adpraut*c00*tem3*(1.-exp(-tem5))
                adtem3 = adtem3+adpraut*c00*amaxcmr*(1.-exp(-tem5))
                adtem5 = adtem5+adpraut*c00*tem3*amaxcmr*exp(-tem5)
                adpraut = 0.
                if (term6 .lt. 50.) then
                  adterm6 = adterm6+adtem5
                  adtem5 = 0.
                else
                  adtem5 = 0.
                endif
                adtem4 = adtem4+2*adterm6*tem4
                adterm6 = 0.
                adamaxcmr = adamaxcmr+adtem4*(cmr*tem3/term5)
                adtem3 = adtem3+adtem4*(amaxcmr*cmr/term5)
                adterm5 = adterm5-adtem4*(amaxcmr*cmr*tem3/(term5*term5))
                adtem4 = 0.
                if (ccr .gt. 0.01) then
                  adccr = adccr+adterm5
                  adterm5 = 0.
                else
                  adterm5 = 0.
                endif

!!old - code below can result in divide by zero
!!                adtem1 = adtem1+adtem3*c1*1./(2.*sqrt(tem1*rdt))* &
!!                     rdt*(1+c2*sqrt(tem2))
!!                adtem2 = adtem2+adtem3*(1.+c1*sqrt(tem1*rdt))*c2* &
!!                     (1./(2.*sqrt(tem2)))
                if (abs(tem1).gt.1.e-20) then 
                   adtem1 = adtem1+adtem3*c1*1./ &
                        (2.*sqrt(tem1*rdt))*rdt*(1+c2*sqrt(tem2))
                else
                   adtem1 = adtem1+adtem3*c1*1./ &
                        (2.*sqrt(1.e-20*rdt))*rdt*(1+c2*sqrt(tem2))
                endif
                if (abs(tem2).gt.1.e-20) then
                   adtem2 = adtem2+adtem3*(1.+c1*sqrt(tem1*rdt))* &
                        c2*(1./(2.*sqrt(tem2)))
                else
                   adtem2 = adtem2+adtem3*(1.+c1*sqrt(tem1*rdt))* &
                        c2*(1./(2.*sqrt(1.e-20)))
                endif
!!new

                adtem3 = 0.
                adterm4 = adterm4+adtem2
                adtem2 = 0.
                if (term3r .lt. 20.) then
                  adterm3r = adterm3r+adterm4
                  adterm4 = 0.
                else
                  adterm4 = 0.
                endif
                if (term2r .gt. 0.) then
                  adterm2r = adterm2r+adterm3r
                  adterm3r = 0.
                else
                  adterm3r = 0.
                endif
                adtin = adtin-adterm2r
                adterm2r = 0.
                adprecrl_0 = adprecrl_0+adtem1
                adprecsl_0 = adprecsl_0+adtem1
                adtem1 = 0.
                adcwmk = adcwmk+adamaxcmr
                adamaxcmr = 0
              endif
              if (ww0 .gt. 0.) then
                adww0 = adww0+adcwmk
                adcwmk = 0.
              else
                adcwmk = 0.
              endif
              adcwmin = adcwmin+adww0
              adww0 = 0.
            else
              adprecsl_0 = adprecsl_0+adprecsl_1
              adprecsl_1 = 0.
              adprecrl_0 = adprecrl_0+adprecrl_1
              adprecrl_1 = 0.
              adcwmin = adcwmin+adww2
              adww2 = 0.
            endif
            if (rq .lt. rhc(k)) then
              adccr = 0.
            else if (rq .ge. us) then
              adccr = 0.
            else
              adrqkll = adrqkll+adccr*(1./(2.*sqrt((us-rqkll)/ &
                   (us-rhc(k))))/(us-rhc(k)))
              adccr = 0.
              if (rq .lt. us) then
                adrq = adrq+adrqkll
                adrqkll = 0.
              else
                adrqkll = 0.
              endif
            endif
            if (qs .le. 1.e-10) then
              adrq = 0.
            else
              adqk = adqk+adrq/qs
              adqs = adqs-adrq*(qk/(qs*qs))
              adrq = 0.
            endif
            if (tmt0 .lt. (-15.)) then
              adqk = adqk+adfi
              adqs = adqs-adfi*rhc(k)
              adfi = 0.
            endif
            if (qs0 .gt. epsq) then
              adqs0 = adqs0+adqs
              adqs = 0.
            else
              adqs = 0.
            endif
            ades = ades+adqs0*(eps/(pres+epsm1*es)-eps*es*epsm1/ &
                 ((pres+epsm1*es)*(pres+epsm1*es)))
            adqs0 = 0.
            call adfpvsx( tin,es,adtin,ades,adjoint )
            ades = 0.
            adtin = adtin+adtmt0
            adtmt0 = 0.
            if (qin .gt. epsq) then
              adqin = adqin+adqk
              adqk = 0.
            else
              adqk = 0.
            endif
          else
            adcwmin = adcwmin+adww2
            adww2 = 0.
            adprecsl_0 = adprecsl_0+adprecsl_3
            adprecsl_3 = 0.
            adprecrl_0 = adprecrl_0+adprecrl_3
            adprecrl_3 = 0.
            adqin = adqin+adq2
            adq2 = 0.
            adtin = adtin+adt2
            adt2 = 0.
          endif
          if (cwmin .gt. climit) then
            adcwmin = adcwmin+adwwn
            adwwn = 0.
          else
            adwwn = 0.
          endif
          if (precsl_0 .gt. 0.) then
            adprecsl_0 = adprecsl_0+adprecsk0
            adprecsk0 = 0.
          else
            adprecsk0 = 0.
          endif
          if (precrl_0 .gt. 0.) then
            adprecrl_0 = adprecrl_0+adprecrk0
            adprecrk0 = 0.
          else
            adprecrk0 = 0.
          endif
          adcwm_in(k,i) = adcwm_in(k,i)+adcwmin
          adcwmin = 0.
          adq_in(k,i) = adq_in(k,i)+adqin
          adqin = 0.
          adt_in(k,i) = adt_in(k,i)+adtin
          adtin = 0.
          adprecsl1k = adprecsl1k+adprecsl_0
          adprecsl_0 = 0.
          adprecrl1k = adprecrl1k+adprecrl_0
          adprecrl_0 = 0.
        end do

        else
          adrn_out(i) = 0.
          do k = km, 1, -1
            if (cwm_in(k,i) .lt. 0.) then
              adcwm_out(k,i) = 0.
              adcwm_in(k,i) = adcwm_in(k,i)-adt_out(k,i)*elwv*rcp
              adt_in(k,i) = adt_in(k,i)+adt_out(k,i)
              adt_out(k,i) = 0.
              adcwm_in(k,i) = adcwm_in(k,i)+adq_out(k,i)
              adq_in(k,i) = adq_in(k,i)+adq_out(k,i)
              adq_out(k,i) = 0.
            else
              adcwm_in(k,i) = adcwm_in(k,i)+adcwm_out(k,i)
              adcwm_out(k,i) = 0.
              adt_in(k,i) = adt_in(k,i)+adt_out(k,i)
              adt_out(k,i) = 0.
              adq_in(k,i) = adq_in(k,i)+adq_out(k,i)
              adq_out(k,i) = 0.
            endif
          end do
        endif
        adrn_out(i) = 0.
        adprecsl1k = 0.
        adprecrl1k = 0.
      end do


      end