*=*=*=*= newcondens.html =*=*=*=*
SUBROUTINE newcondens

SUBROUTINE newcondens


      SUBROUTINE newcondens(ngrid,nlayer,ptimestep,
     $                  pcapcal,pplay,pplev,ptsrf,pt,
     $                  pphi,pdt,pdu,pdv,pdtsrf,pu,pv,paerosol,
     $                  piceco2,psolaralb,pemisurf,
     $                  pdtc,pdtsrfc,pdpsrf,pduc,pdvc,
     $                  pdaerosol)

       IMPLICIT NONE
c=======================================================================
c   subject:
c   --------
c   Condensation/sublimation of CO2 ice on the ground and in the
c   atmosphere
c   (Scheme described in Forget et al., Icarus, 1998)
c
c   author:   Francois Forget     1994-1996
c   ------
c
c   input:
c   -----
c    ngrid                 nombre de points de verticales
c                          (toutes les boucles de la physique sont au
c                          moins vectorisees sur ngrid)
c    nlayer                nombre de couches
c    pplay(ngrid,nlayer)   Pressure levels
c    pplev(ngrid,nlayer+1) Pressure levels
c    pt(ngrid,nlayer)      temperature (en K)
c    ptsrf(ngrid)          temperature de surface
c
c                    \
c    pdt(ngrid,nlayermx)   \  derivee temporelle physique avant condensation
c                     /  ou sublimation pour  pt,ptsrf
c    pdtsrf(ngrid)   /
c
c   output:
c   -------
c
c    pdpsrf(ngrid)   \  derivee temporelle physique (contribution de
c    pdtc(ngrid,nlayermx) /  la condensation ou sublimation) pour Ps,pt,ptsrf
c    pdtsrfc(ngrid) /
c
c   Entree/sortie
c   -------------
c
c    piceco2(ngrid) :      quantite de glace co2 au sol (kg/m2)
c    psolaralb(ngrid,2) :  albedo au sol
c    pemisurf(ngrid)     :  emissivite du sol

c
c=======================================================================
c
c    0.  Declarations :
c    ------------------
c
#include "dimensions.h"
#include "dimphys.h"
#include "comcstfi.h"
#include "surfdat.h"
#include "comgeomfi.h"
#include "comvert.h"


c-----------------------------------------------------------------------
c    Arguments :
c    ---------
      INTEGER ngrid,nlayer

      REAL ptimestep
      REAL pplay(ngrid,nlayer),pplev(ngrid,nlayer+1)
      REAL pcapcal(ngrid)
      REAL pt(ngrid,nlayer)
      REAL ptsrf(ngrid)
      REAL pphi(ngrid,nlayer)
      REAL pdt(ngrid,nlayer),pdtsrf(ngrid),pdtc(ngrid,nlayer)
      REAL pdtsrfc(ngrid),pdpsrf(ngrid)
      REAL piceco2(ngrid),psolaralb(ngrid,2),pemisurf(ngrid)

      REAL pu(ngrid,nlayer) ,  pv(ngrid,nlayer)
      REAL pdu(ngrid,nlayer) , pdv(ngrid,nlayer)
      REAL pduc(ngrid,nlayer) , pdvc(ngrid,nlayer)
      REAL paerosol(ngridmx,nlayermx,5)
      REAL pdaerosol(ngridmx,nlayermx)
c
c    Local variables :
c    -----------------
      INTEGER l,ig,icap

c
      REAL zt(ngridmx,nlayermx)
      REAL zcpi
      REAL ztcond (ngridmx,nlayermx)
      REAL ztcondsol(ngridmx)
      REAL zdiceco2(ngridmx)
      REAL zcondicea(ngridmx,nlayermx), zcondices(ngridmx)
      REAL zfallice(ngridmx,nlayermx+1) , zfallheat
      REAL zmflux(nlayermx+1)
      REAL zu(nlayermx),zv(nlayermx)
      REAL ztsrf(ngridmx)
      REAL ztc(nlayermx), ztm(nlayermx+1)
      REAL zum(nlayermx+1) , zvm(nlayermx+1)
      REAL masse(nlayermx),w(nlayermx+1)
      LOGICAL condsub(ngridmx)

c variable speciale diagnostique
      real tconda1(ngridmx,nlayermx)
      real tconda2(ngridmx,nlayermx)
c     REAL zdiceco2a(ngridmx) ! for diagnostic only
      real zdtsig (ngridmx,nlayermx)
      real zdt (ngridmx,nlayermx)



c   local saved variables
      REAL emisref(ngridmx)
      REAL latcond,tcond1mb
      REAL acond,bcond,ccond,cpice
      SAVE emisref
      SAVE latcond,acond,bcond,ccond,cpice


      LOGICAL firstcall
      SAVE firstcall
      REAL SSUM
      EXTERNAL SSUM

      common/scratch/zt,ztcond,zcondicea,zfallice,tconda1
     ,    ,tconda2,zdtsig,zdt,zu,zv


      DATA latcond,tcond1mb/5.9e5,136.27/
      DATA cpice /1000./
      DATA firstcall/.true./

c----------------------------------------------------------------------

c   Initialisation
c   --------------
c
      IF (firstcall) THEN
         bcond=1./tcond1mb
         ccond=cpp/(g*latcond)
         acond=r/latcond
         firstcall=.false.
         PRINT*,'In newcondens:Tcond(P=1mb)=',tcond1mb,' Lcond=',latcond
         PRINT*,'acond,bcond,ccond',acond,bcond,acond
      ENDIF
      zcpi=1./cpp
c
c======================================================================
c    Calcul of CO2 condensation sublimation
c    ============================================================
c
c    Used variable :
c       piceco2(ngrid)   :  amount of co2 ice on the ground (kg/m2)
c       zcondicea(ngrid,l):  condensation rate in layer  l  (kg/m2/s)
c       zcondices(ngrid):  condensation rate on the ground (kg/m2/s)
c       zfallice(ngrid,l):amount of ice falling from layer l (kg/m2/s)
c
c       pdtc(ngrid,nlayermx)  : dT/dt due to cond/sub
c
c
c     Tendancies set to 0 (except pdtc)
c     -------------------------------------
      DO l=1,nlayer
         DO ig=1,ngrid
           zcondicea(ig,l) = 0.
           zfallice(ig,l) = 0.
           pduc(ig,l)  = 0
           pdvc(ig,l)  = 0
           pdaerosol(ig,l) = 0
         END DO
      END DO

      DO ig=1,ngrid
         zfallice(ig,nlayer+1) = 0.
         zcondices(ig) = 0.
         pdtsrfc(ig) = 0.
         pdpsrf(ig) = 0.
         condsub(ig) = .false.
         zdiceco2(ig) = 0.
      ENDDO
      zfallheat=0

c     *************************
c     ATMOSPHERIC CONDENSATION
c     *************************

c     forecast of atmospheric temperature zt and frost temperature ztcond
c     --------------------------------------------------------------------

      DO l=1,nlayer
         DO ig=1,ngrid
            zt(ig,l)=pt(ig,l)+ pdt(ig,l)*ptimestep
            ztcond(ig,l)=1./(bcond-acond*log(.0095*pplay(ig,l)))
         ENDDO
      ENDDO


c      Condensation/sublimation in the atmosphere
c      ------------------------------------------
c      (calcul of zcondicea , zfallice and pdtc)
c
      DO l=nlayer , 1, -1
         DO ig=1,ngrid
           pdtc(ig,l)=0.
           IF((zt(ig,l).LT.ztcond(ig,l)).or.(zfallice(ig,l+1).gt.0))THEN
               condsub(ig)=.true.
               IF (zfallice(ig,l+1).gt.0) then
                 zfallheat=zfallice(ig,l+1)*
     &           (pphi(ig,l+1)-pphi(ig,l) +
     &          cpice*(ztcond(ig,l+1)-ztcond(ig,l)))/latcond
               ELSE
                   zfallheat=0.
               ENDIF
               pdtc(ig,l)=(ztcond(ig,l) - zt(ig,l))/ptimestep
               zcondicea(ig,l)=(pplev(ig,l)-pplev(ig,l+1))
     &                        *ccond*pdtc(ig,l)- zfallheat
c              Case when the ice from above sublimes entirely
c              """""""""""""""""""""""""""""""""""""""""""""""
               IF (zfallice(ig,l+1).lt.- zcondicea(ig,l)) then
                  pdtc(ig,l)=(-zfallice(ig,l+1)+zfallheat)/
     &              (ccond*(pplev(ig,l)-pplev(ig,l+1)))
                  zcondicea(ig,l)= -zfallice(ig,l+1)
               END IF

               zfallice(ig,l) = zcondicea(ig,l)+zfallice(ig,l+1)
            END IF
         ENDDO
      ENDDO

c     *************************
c       SURFACE  CONDENSATION
c     *************************

c     forecast of ground temperature ztsrf and frost temperature ztcondsol
c     --------------------------------------------------------------------
      DO ig=1,ngrid
         ztcondsol(ig)=1./(bcond-acond*log(.0095*pplev(ig,1)))
         ztsrf(ig) = ptsrf(ig) + pdtsrf(ig)*ptimestep
      ENDDO

c
c      Condensation/sublimation on the ground
c      --------------------------------------
c      (calcul of zcondices , pdtsrfc)
c
      DO ig=1,ngrid
         IF(ig.GT.ngrid/2+1) THEN
            icap=2
         ELSE
            icap=1
         ENDIF

c
c        Loop on where we have  condensation/ sublimation
         IF ((ztsrf(ig) .LT. ztcondsol(ig)) .OR.   ! ground cond
     $       (zfallice(ig,1).NE.0.) .OR.            ! falling snow
     $      ((ztsrf(ig) .GT. ztcondsol(ig)) .AND.  ! ground sublim.
     $      ((piceco2(ig)+zfallice(ig,1)*ptimestep) .NE. 0.))) THEN
            condsub(ig) = .true.

            IF (zfallice(ig,1).gt.0) then
                 zfallheat=zfallice(ig,1)*
     &           (pphi(ig,1)- phisfi(ig) +
     &           cpice*(ztcond(ig,1)-ztcondsol(ig)))/(latcond*ptimestep)
            ELSE
                 zfallheat=0.
            ENDIF

c           condensation or partial sublimation of CO2 ice
c           """""""""""""""""""""""""""""""""""""""""""""""
            zcondices(ig)=pcapcal(ig)*(ztcondsol(ig)-ztsrf(ig))
     &      /(latcond*ptimestep)         - zfallheat
            pdtsrfc(ig) = (ztcondsol(ig) - ztsrf(ig))/ptimestep

c           If the entire CO_2 ice layer sublime
c           """"""""""""""""""""""""""""""""""""""""""""""""""""
c           (including what has just condensed in the atmosphere)

            IF((piceco2(ig)/ptimestep+zfallice(ig,1)).LE.
     &          -zcondices(ig))THEN
               zcondices(ig) = -piceco2(ig)/ptimestep - zfallice(ig,1)
               pdtsrfc(ig)=(latcond/pcapcal(ig))*
     &         (zcondices(ig)+zfallheat)
            END IF

c           Changing CO2 ice amount and pressure :
c           """"""""""""""""""""""""""""""""""""

            zdiceco2(ig) = zcondices(ig) + zfallice(ig,1)
            piceco2(ig) = piceco2(ig) + zdiceco2(ig)*ptimestep
            pdpsrf(ig) = -zdiceco2(ig)*g

            IF(ABS(pdpsrf(ig)*ptimestep).GT.pplev(ig,1)) THEN
               PRINT*,'STOP in condens'
               PRINT*,'condensing more than total mass'
               PRINT*,'Grid point ',ig
               PRINT*,'Ps = ',pplev(ig,1)
               PRINT*,'d Ps = ',pdpsrf(ig)
               STOP
            ENDIF
         END IF
      ENDDO

c  Surface albedo and emissivity of the ground below the snow (emisref)
c  --------------------------------------------------------------------
         do ig =1,ngrid
           IF(ig.GT.ngrid/2+1) THEN
              icap=2
           ELSE
              icap=1
           ENDIF

           if(.not.piceco2(ig).ge.0.) THEN
              print*,'WARNING newcondens piceco2(',ig,')=',
     .          piceco2(ig)
               piceco2(ig)=0.
           endif
           if (piceco2(ig).gt.0) then
               psolaralb(ig,1) = albedice(icap)
               psolaralb(ig,2) = albedice(icap)
               emisref(ig) = emisice(icap)
           else
               psolaralb(ig,1) = albedodat(ig)
               psolaralb(ig,2) = albedodat(ig)
               emisref(ig) = emissiv
               pemisurf(ig) = emissiv
           end if
         end do


c ***************************************************************
c  Correction to account for redistribution between sigma layers
c  when changing surface pressure (and warming/cooling
c  of the CO2 currently changing phase).
c *************************************************************

      DO ig=1,ngrid
        if (condsub(ig)) then
           do l=1,nlayer
             ztc(l) = zt(ig,l) + pdtc(ig,l)*ptimestep
             zu(l)=pu(ig,l)+ pdu(ig,l)*ptimestep
             zv(l)=pv(ig,l)+ pdv(ig,l)*ptimestep
           end do

c  Mass fluxes through the sigma levels (kg.m-2.s-1)  (>0 when up)
c  """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
           zmflux(1) = -zcondices(ig)
           DO l=1,nlayer
              zmflux(l+1) = zmflux(l)
     &        -zcondicea(ig,l) + dsig(l)*(zfallice(ig,1)-zmflux(1))
           END DO

c Mass of each layer
c ------------------
           DO l=1,nlayer
              masse(l)=(pplev(ig,l) - pplev(ig,l+1))/g
           END DO


c  Corresponding fluxes for T,U,V,Q
c  """"""""""""""""""""""""""""""""

c           averaging operator for TRANSPORT
c           """"""""""""""""""""""""""""""""
            ztm(1)= ztsrf(ig) + pdtsrfc(ig)*ptimestep
            zum(1)= 0
            zvm(1)= 0

c           Centered scheme (no more used):
c           DO l=2,nlayer
c               ztm(l)= 0.5*(ztc(l-1)+ztc(l))
c               zum(l)= 0.5*(zu(l-1)+zu(l))
c               zvm(l)= 0.5*(zv(l-1)+zv(l))
c            END DO

c           Van Leer scheme:
            DO l=1,nlayer+1
                w(l)=-zmflux(l)*ptimestep
            END DO
            call vl1d(ztc,2.,masse,w,ztm)
            call vl1d(zu,2.,masse,w,zum)
            call vl1d(zv,2.,masse,w,zvm)

c           Surface condensation affects low winds
            if (zmflux(1).lt.0) then
                zum(1)= zu(1) *  (w(1)/masse(1))
                zvm(1)= zv(1) *  (w(1)/masse(1))
                if (w(1).gt.masse(1)) then ! ensure numerical stability
                  zum(1)= (zu(1)-zum(2))*masse(1)/w(1) + zum(2)
                  zvm(1)= (zv(1)-zvm(2))*masse(1)/w(1) + zvm(2)
                end if
            end if

            ztm(nlayer+1)= ztc(nlayer) ! should not be used, but...
            zum(nlayer+1)= zu(nlayer) ! should not be used, but...
            zvm(nlayer+1)= zv(nlayer) ! should not be used, but...

c           Tendancies on T, U, V
c           """""""""""""""""""""
            DO l=1,nlayer

c            Tendancies on T
              zdtsig(ig,l)  =  (1/masse(l)) *
     &      ( zmflux(l)*(ztm(l) - ztc(l))
     &      - zmflux(l+1)*(ztm(l+1) - ztc(l))
     &      + zcondicea(ig,l)*(ztcond(ig,l)-ztc(l))  )
              pdtc(ig,l) =  pdtc(ig,l) + zdtsig(ig,l)

c           Tendancies on U
              pduc(ig,l) =   (1/masse(l)) *
     &      ( zmflux(l)*(zum(l) - zu(l))
     &      - zmflux(l+1)*(zum(l+1) - zu(l)) )


c           Tendancies on V
              pdvc(ig,l) =   (1/masse(l)) *
     &      ( zmflux(l)*(zvm(l) - zv(l))
     &      - zmflux(l+1)*(zvm(l+1) - zv(l)) )

            END DO
        END IF
      END DO


c ***************************************************************
c   CO2 snow / clouds scheme
c ***************************************************************

      call co2snow(ngrid,nlayer,ptimestep,emisref,condsub,pplev,
     &        paerosol,zcondicea,zcondices,zfallice,pemisurf,
     &             pdaerosol)

c ***************************************************************
c Ecriture des diagnostiques
c ***************************************************************
c     DO l=1,nlayer
c        DO ig=1,ngrid
c         Taux de cond en kg.m-2.pa-1.s-1
c          tconda1(ig,l)=zcondicea(ig,l)/(pplev(ig,l)-pplev(ig,l+1))
c          Taux de cond en kg.m-3.s-1
c          tconda2(ig,l)=tconda1(ig,l)*pplay(ig,l)*g/(r*pt(ig,l))
c        END DO
c     END DO
c     DO ig=1,ngrid
c          zdiceco2a(ig)=zfallice(ig,1)
c     END DO
c     call WRITEDRSFI(ngridmx,44,'diceco2a',
c    &'Taux de cond atm','g.m-2.s-1',2,zdiceco2a)
c     call WRITEDRSFI(ngridmx,44,'diceco2s',
c    &'Taux de cond sol','g.m-2.s-1',2,zcondices)
c     call WRITEDRSFI(ngridmx,44,'tconda1',
c    &'Taux de condensation CO2 atmospherique /Pa',
c    & 'kg.m-2.Pa-1.s-1',3,tconda1)
c     call WRITEDRSFI(ngridmx,44,'tconda2',
c    &'Taux de condensation CO2 atmospherique /m',
c    & 'kg.m-3.s-1',3,tconda2)
c     call WRITEDRSFI(ngridmx,44,'fallice',
c    &'falling ice', 'kg.m-2.s-1',3,zfallice)
c     call WRITEDRSFI(ngridmx,44,'emisref',
c    &'emisref',' ',2,emisref)

c     call WRITEDRSFI(ngridmx,44,'pdtsrfc',
c    &'Tendance Tsol','K.s-1',2, pdtsrfc)
c     call WRITEDRSFI(ngridmx,44,'pdtc',
c    &'Tendance T due a cond', 'K.s-1',3,pdtc)
c     call WRITEDRSFI(ngridmx,44,'zdtsig',
c    &'Tendance T due a sigma aj.', 'K.s-1',3,zdtsig)
c     call WRITEDRSFI(ngridmx,44,'ducond ',
c    &'Tendance U due a cond', 'm.s-2',3,pduc)
c     call WRITEDRSFI(ngridmx,44,'dvcond',
c    &'Tendance V due a cond', 'm.s-2',3,pdvc)

      return
      end
c ***************************************************************** *=*=*=*= vl1d.html =*=*=*=*
SUBROUTINE vl1d

SUBROUTINE vl1d


      SUBROUTINE vl1d(q,pente_max,masse,w,qm)
c
c
c     Operateur de moyenne inter-couche pour calcul de transport type
c     Van-Leer " pseudo amont " dans la verticale
c    q,w sont des arguments d'entree  pour le s-pg ....
c    masse : masse de la couche Dp/g
c    w : masse d'atm ``transferee'' a chaque pas de temps (kg.m-2)
c    pente_max = 2 conseillee
c
c
c   --------------------------------------------------------------------
      IMPLICIT NONE

#include "dimensions.h"

c
c
c
c   Arguments:
c   ----------
      real masse(llm),pente_max
      REAL q(llm),qm(llm+1)
      REAL w(llm+1)
c
c      Local
c   ---------
c
      INTEGER l
c
      real dzq(llm),dzqw(llm),adzqw(llm),dzqmax
      real sigw, Mtot, MQtot
      integer m
c     integer ismax,ismin


c    On oriente tout dans le sens de la pression
c     W > 0 WHEN DOWN !!!!!!!!!!!!!

      do l=2,llm
            dzqw(l)=q(l-1)-q(l)
            adzqw(l)=abs(dzqw(l))
      enddo

      do l=2,llm-1
            if(dzqw(l)*dzqw(l+1).gt.0.) then
                dzq(l)=0.5*(dzqw(l)+dzqw(l+1))
            else
                dzq(l)=0.
            endif
            dzqmax=pente_max*min(adzqw(l),adzqw(l+1))
            dzq(l)=sign(min(abs(dzq(l)),dzqmax),dzq(l))
      enddo

         dzq(1)=0.
         dzq(llm)=0.

       do l = 1,llm-1

c         Regular scheme (transfered mass < layer mass)
c         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
          if(w(l+1).gt.0. .and. w(l+1).le.masse(l+1)) then
             sigw=w(l+1)/masse(l+1)
             qm(l+1)=(q(l+1)+0.5*(1.-sigw)*dzq(l+1))
          else if(w(l+1).le.0. .and. -w(l+1).le.masse(l)) then
             sigw=w(l+1)/masse(l)
             qm(l+1)=(q(l)-0.5*(1.+sigw)*dzq(l))

c         Extended scheme (transfered mass > layer mass)
c         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
          else if(w(l+1).gt.0.) then
             m=l+1
             Mtot = masse(m)
             MQtot = masse(m)*q(m)
             do while ((m.lt.llm).and.(w(l+1).gt.(Mtot+masse(m+1))))
                m=m+1
                Mtot = Mtot + masse(m)
                MQtot = MQtot + masse(m)*q(m)
             end do
             if (m.lt.llm) then
                sigw=(w(l+1)-Mtot)/masse(m+1)
                qm(l+1)= (1/w(l+1))*(MQtot + (w(l+1)-Mtot)*
     &          (q(m+1)+0.5*(1.-sigw)*dzq(m+1)) )
             else
                w(l+1) = Mtot
                qm(l+1) = Mqtot / Mtot
                write(*,*) 'top layer is disapearing !'
                stop
             end if
          else      ! if(w(l+1).lt.0)
             m = l-1
             Mtot = masse(m+1)
             MQtot = masse(m+1)*q(m+1)
             do while ((m.gt.0).and.(-w(l+1).gt.(Mtot+masse(m))))
                m=m-1
                Mtot = Mtot + masse(m+1)
                MQtot = MQtot + masse(m+1)*q(m+1)
             end do
             if (m.gt.0) then
                sigw=(w(l+1)+Mtot)/masse(m)
                qm(l+1)= (-1/w(l+1))*(MQtot + (-w(l+1)-Mtot)*
     &          (q(m)-0.5*(1.+sigw)*dzq(m))  )
             else
                qm(l+1)= (-1/w(l+1))*(MQtot + (-w(l+1)-Mtot)*qm(1))
             end if
          end if
       enddo

c     boundary conditions (not used in newcondens !!)
c         qm(llm+1)=0.
c         if(w(1).gt.0.) then
c            qm(1)=q(1)
c         else
c           qm(1)=0.
c         end if

      return
      end