MODULE my_smom_mod 1

  implicit none

  private
  public :: mysmom_main

  contains

!_______________________________________________________________________________________!

  SUBROUTINE MYSMOM_MAIN(W_omega,T,Q,QC,QR,QI,QN,QG,QH,PS,TM,QM,QCM,QRM,QIM,QNM,QGM,    & 1,59
     QHM,PSM,S,RT_rn1,RT_rn2,RT_fr1,RT_fr2,RT_sn1,RT_sn2,RT_sn3,RT_pe1,RT_pe2,RT_peL,   &
     RT_snd,GZ,T_TEND,Q_TEND,QCTEND,QRTEND,QITEND,QNTEND,QGTEND,QHTEND,dt,NI,N,NK,J,    &
     KOUNT,CCNtype,dzMIN,Dm_c,Dm_r,Dm_i,Dm_s,Dm_g,Dm_h,ZET,ZEC,SLW,VIS,VIS1,VIS2,VIS3,  &
     h_CB,h_ML1,h_ML2,h_SN)

  use my_fncs_mod
  use my_sedi_mod

  IMPLICIT NONE

!CALLING PARAMETERS:
  integer,                intent(in)    :: NI,NK,N,J,KOUNT,CCNtype
  real,                   intent(in)    :: dt,dzMIN
  real, dimension(:),     intent(in)    :: PS,PSM
  real, dimension(:),     intent(out)   :: h_CB,h_ML1,h_ML2,h_SN
  real, dimension(:),     intent(out)   :: RT_rn1,RT_rn2,RT_fr1,RT_fr2,RT_sn1,RT_sn2, &
                                           RT_sn3,RT_pe1,RT_pe2,RT_peL,RT_snd,ZEC
  real, dimension(:,:),   intent(in)    :: W_omega,S,GZ
  real, dimension(:,:),   intent(inout) :: T,Q,QC,QR,QI,QN,QG,QH,TM,QM,QCM,QRM,QIM,   &
                                           QNM,QGM,QHM
  real, dimension(:,:),   intent(out)   :: T_TEND,QCTEND,QRTEND,QITEND,QNTEND,QGTEND, &
                                           QHTEND,Q_TEND,ZET,Dm_c,Dm_r,Dm_i,Dm_s,     &
                                           Dm_g,Dm_h,SLW,VIS,VIS1,VIS2,VIS3

!_______________________________________________________________________________________
!                                                                                       !
!                    Milbrandt-Yau Multimoment Bulk Microphysics Scheme                 !
!                       - six-category, single-moment version --                        !
!_______________________________________________________________________________________!
!
!  Author:
!       J. Milbrandt, McGill University (August 2004)
!
!  Revision:
!
!  001  J. Milbrandt  (Dec 2006) - Converted the full Milbrandt-Yau (2005) multimoment
!         (RPN)                    scheme to an efficient single-moment version
!  002  J. Milbrandt  (Dec 2007) - Modifications for implemntation into GEM-LAM-2.5
!                                   + new precipitation types (RT_x),
!                                   + unified call to sedimentation for all categories
!                                   + bug corrections
!  003  J. Milbrandt  (Dec 2007) - Added explicit interface for call in vkuocon6.ftn.
!
!
!  Object:
!          Computes changes to the temperature, water vapor mixing ratio, and the
!          mixing ratios of six hydrometeor species resulting from cloud microphysical
!          interactions at saturated grid points. Surface precipitation rates from each
!          sedimenting hydrometeor categories are also computed.
!
!          This subroutine and the associated modules form the single-moment version of
!          the multimoment bulk microphysics package, described in the references below.
!
!  References:   Milbrandt and Yau, (2005a), J. Atmos. Sci., vol.62, 3051-3064
!                --------- and ---, (2005b), J. Atmos. Sci., vol.62, 3065-3081
!                (and references therein)
!
!  Please report bugs to:  jason.milbrandt@ec.gc.ca
!_______________________________________________________________________________________
!  Package version:   2.12.2                                                            !
!  Last modified  :   2009-01-20                                                        !
!_______________________________________________________________________________________!
!
! Arguments:         Description:                                         Units:
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -!
!            - Input -
!
! NI                 number of x-dir points (local subdomain)
! NK                 number of vertical levels
! N
! J                  y-dir index (local subdomain)
! KOUNT              current model time step number
! dt                 model time step                                      [s]
! CCNtype            switch for airmass type
!                      1 = maritime                   --> N_c =  80 cm-3
!                      2 = continental 1              --> N_c = 200 cm-3
!                      3 = continental 2  (polluted)  --> N_c = 500 cm-3
!                      4 = land-sea-mask-dependent (TBA)
! dzMIN              thickness of lowest level, used for sedimentation    [m]
! W_omega            vertical velocity                                    [Pa s-1]
! S                  sigma
! GZ                 geopotential height                                  [m]
!
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -!
!            - Input/Output -
!
! T                  air temperature at time (t*)                         [K]
! TM                 air temperature at time (t-dt)                       [K]
! Q                  water vapor mixing ratio at (t*)                     [kg kg-1]
! QM                 water vapor mixing ratio at (t-dt)                   [kg kg-1]
! PS                 surface pressure at time (t*)                        [Pa]
! PSM                surface pressure at time (t-dt)                      [Pa]
!
!  For x = (C,R,I,N,G,H):  C = cloud
!                          R = rain
!                          I = ice (pristine) [except 'NY', not 'NI']
!                          N = snow
!                          G = graupel
!                          H = hail
!
! Q(x)               mixing ratio for hydrometeor x at (t*)               [kg kg-1]
! Q(x)M              mixing ratio for hydrometeor x at (t-dt)             [kg kg-1]
! N(x)               total number concentration for hydrometeor x  (t*)   [m-3]
! N(x)M              total number concentration for hydrometeor x  (t-dt) [m-3]
!
! Note:  The arrays "VM" (e.g. variables TM,QM,QCM etc.) are declared as INTENT(INOUT)
!        such that their values are modified in the code [VM = 0.5*(VM + V)].
!        This is to approxiate the values at time level (t), which are needed by
!        this routine but are unavailable to the PHYSICS.  The new values are discared
!        by the calling routine ('vkuocon6.ftn').  However, care should be taken with
!        interfacing with other modelling systems.  For GEM/MC2, it does not matter if
!        VM is modified since the calling module passes back only the tendencies
!        (VTEND) to the model.

!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -!
!            - Output -
!
! Q_TEND             tendency for water vapor mixing ratio                [kg kg-1 s-1]
! T_TEND             tendency for air temperature                         [K s-1]
! Q(x)TEND           tendency for mixing ratio for hydrometeor x          [kg kg-1 s-1]
! N(x)TEND           tendency for number concentration for hydrometeor x  [m-3 s-1]
! Dm_(x)             mean-mass diameter for hydrometeor x                 [m]
! RT_rn1             precipitation rate (at sfc) of liquid rain           [m+3 m-2 s-1]
! RT_rn2             precipitation rate (at sfc) of liquid drizzle        [m+3 m-2 s-1]
! RT_fr1             precipitation rate (at sfc) of freezing rain         [m+3 m-2 s-1]
! RT_fr2             precipitation rate (at sfc) of freezing drizzle      [m+3 m-2 s-1]
! RT_sn1             precipitation rate (at sfc) of ice crystals (liq-eq) [m+3 m-2 s-1]
! RT_sn2             precipitation rate (at sfc) of snow    (liq-equiv)   [m+3 m-2 s-1]
! RT_sn3             precipitation rate (at sfc) of graupel (liq-equiv)   [m+3 m-2 s-1]
! RT_snd             precipitation rate (at sfc) of snow    (frozen)      [m+3 m-2 s-1]
! RT_pe1             precipitation rate (at sfc) of ice pellets (liq-eq)  [m+3 m-2 s-1]
! RT_pe2             precipitation rate (at sfc) of hail (total; liq-eq)  [m+3 m-2 s-1]
! RT_peL             precipitation rate (at sfc) of hail (large only)     [m+3 m-2 s-1]
! SLW                supercooled liquid water content                     [kg m-3]
! VIS                visibility resulting from fog, rain, snow            [m]
! VIS1               visibility component through liquid cloud (fog)      [m]
! VIS2               visibility component through rain                    [m]
! VIS3               visibility component through snow                    [m]
! ZET                total equivalent radar reflectivity                  [dBZ]
! ZEC                composite (column-max) of ZET                        [dBZ]
!_______________________________________________________________________________________!


! LOCAL VARIABLES:

 !Variables to count active grid points:
  logical :: log1,log2,log3,log4,doneK,rainPresent,activePoint(ni,nk),calcDiag,     &
             CB_found,ML_found,SN_found
  integer :: i,k,niter,ll,start,ktop_sedi

  real    :: tmp1,tmp2,tmp3,tmp4,tmp5,tmp6,tmp7,tmp8,tmp9,tmp10,                    &
       VDmax,NNUmax,X,D,DEL,QREVP,ES,NuDEPSOR,NuCONTA,NuCONTB,NuCONTC,iMUkin,       &
       NuCONT,GG,Na,Tcc,F1,F2,Kdiff,PSIa,Kn,source,sink,sour,ratio,qvs0,            &
       DELqvs,ft,esi,Si,Simax,Vq,Vn,Vz,LAMr,No_i,iEih,iEsh,N_r,N_i,N_s,N_g,N_h,     &
       iABi,ABw,VENTr,VENTs,VENTg,VENTi,VENTh,Cdiff,Ka,MUdyn,MUkin,DEo,             &
       gam,ScTHRD,Tc,mi,ri,ff,Ec,Ntr,Dho,DMrain,Ech,DMice,DMsnow,DMgrpl,DMhail,     &
       ssat,Swmax,dey,Esh,Eii,Eis,Ess,Eig,Eih,FRAC,JJ,Dirg,Dirh,Dsrs,Dsrg,Dsrh,     &
       Dgrg,Dgrh,SIGc,L,TAU,DrAUT,DrINIT,Di,Ds,Dg,Dh,qFact,nFact,Ki,Rz,GR37,        &
       vr0,vi0,vs0,vg0,vh0,Dc,Dr,QCLcs,QCLrs,QCLis,QCLcg,QCLrg,QCLig,des_mlt,       &
       QCLch,QCLrh,QCLsh,QMLir,QMLsr,QMLgr,QMLhr,QCLih,QCNig,QVDvg,QVDvh,QSHhr,     &
       QFZci,QNUvi,QCLci,QVDvi,QCNis,QCNis1,QCNis2,QCLir,QCLri,QCNsg,QCLsr,QCNgh,   &
       QCLgr,QHwet,QVDvs,QFZrh,QIMsi,QIMgi,NMLhr,NCNgh,NVDvh,NCLir,NCLri,NCLrh,     &
       NCLch,NCLsr,NCLirg,NCNig,NCLirh,NrFZrh,NhFZrh,NCLsrs,NCLsrg,NCLsrh,NCLgrg,   &
       NCLgrh,NVDvg,NMLgr,NiCNis,NsCNis,NVDvs,NMLsr,NCLsh,NCLss,NNUvi,NVDvi,        &
       NCLis,NCLig,NCLih,NMLir,NCLrs,NCNsg,NCLcs,NCLcg,NCLci,NIMsi,NIMgi,NIMii,     &
       NCLgr,NCLrg,NSHhr,RCAUTR,RCACCR,CCACCR,CCSCOC,CCAUTR,CRSCOR,ALFx,GX2,GX5,    &
       LAMx,iLAMx,iLAMxB0,Dx,ffx,iDE,iGH31,iLAMc,iNCM,iNRM,iNYM,iNNM,iNGM,          &
       iLAMg,iLAMg2,iLAMgB0,iLAMgB1,iLAMgB2,iLAMh,iLAMhB0,iLAMhB1,iLAMhB2,iNHM,     &
       iLAMi,iLAMi2,iLAMi3,iLAMi4,iLAMi5,iLAMiB0,iLAMiB1,iLAMiB2,iLAMr6,iLAMh2,     &
       iLAMs,iLAMs2,iLAMsB0,iLAMsB1,iLAMsB2,iLAMr,iLAMr2,iLAMr3,iLAMr4,iLAMr5,      &
       iLAMc2,iLAMc3,iLAMc4,iLAMc5,iLAMc6,iQCM,iQRM,iQIM,iQNM,iQGM,iQHM,Nt_i,No_s


 !Variables that only need to be calulated on the first step (and saved):
 ! Note: should be declared with SAVE attribute (but not yet tested for openMP)
  real :: idt,iMUc,cmr,cmi,cms,cmg,cmh,icmr,icmi,icms,icmg,icmh,                    &
       GC1,iGC1,GC2,GC3,GC4,GC5,iGC5,GC6,GC7,GC8,GC11,GC12,GC13,GC14,GC15,          &
       GR1,GR3,GR13,GR14,GR15,GR17,GR31,iGR31,GR32,GR33,GR34,GR35,GR36,GI4,         &
       GI6,GI20,GI21,GI22,GI31,GI32,GI33,GI34,GI35,iGI31,GI11,GI36,GI37,GI40,       &
       GS09,GS11,GS12,GS13,iGS20,GS31,iGS31,GS32,GS33,GS34,GS35,GS36,GS40,          &
       GG09,GG11,GG12,GG13,GG31,iGG31,GG32,GG33,GG34,GG35,GG36,GG40,iGG99,          &
       GH09,GH11,GH12,GH13,GH31,GH32,GH33,iGH34,GH40,iGH99,GR50,GG50,GH50,          &
       cexr4,cexr5,cexr6,cexr9,icexr9,ckQr1,ckQr2,ckQr3,ckQi1,ckQi2,ckQi3,ckQi4,    &
       icexi9,ckQs1,ckQs2,ckQs4,cexs1,cexs2,ckQg1,ckQg2,ckQg4,ckQh1,ckQh2,ckQh4,    &
       LCP,LFP,LSP,ck5,ck6,PI2,PIov4,PIov6,CHLS,iCHLF,cxr,cxi,Gzr,Gzi,Gzs,Gzg,Gzh,  &
       icexc9,cexr1,cexr2,cexr3,idew,N_c

!BLG (box-Lagrangian) sedimentation variables:
  integer :: nnn,npassr,npassi,npasss,npassg,npassh,a,counter
  real    :: dtr,dti,dtg,dts,dth,tmp,VqMax,VzMax,cr6,ci6,cg6,cs6,ch6,ck7,CoMax,     &
             CoMAXr,CoMAXi,CoMAXs,CoMAXg,CoMAXh
  logical :: LOCALLIM,slabHASmass

! Size distribution parameters:
  real, parameter :: MUc      =  3.     !shape parameter 1 for cloud
  real, parameter :: alpha_c  =  1.     !shape parameter 2 for cloud
  real, parameter :: alpha_r  =  0.     !shape parameter for rain
  real, parameter :: alpha_i  =  0.     !shape parameter for ice
  real, parameter :: alpha_s  =  0.     !shape parameter for snow
  real, parameter :: alpha_g  =  0.     !shape parameter for graupel
  real, parameter :: alpha_h  =  0.     !shape parameter for hail
  real, parameter :: No_r     =  1.e+7  !intercept parameter for rain    [m-4]
!!real, parameter :: No_s     =  1.e+7  !intercept parameter for snow    [m-4]
!!  (diagnostic No_s=f(T), is used rather than a fixed value)
  real, parameter :: No_g     =  4.e+6  !intercept parameter for graupel [m-4]
  real, parameter :: No_h     =  1.e+5  !intercept parameter for hail    [m-4]
  !------------------------------------!
  ! Symbol convention: (dist. params.) ! MY05: Milbrandt & Yau, 2005a,b (JAS)
  !       MY05    F94       CP00       ! F94:  Ferrier, 1994            (JAS)
  !       ------  --------  ------     ! CP00: Cohard & Pinty, 2000a,b  (QJGR)
  !       ALFx    ALPHAx    MUx-1      !
  !       MUx     (1)       ALPHAx     !
  !       ALFx+1  ALPHAx+1  MUx        !
  !------------------------------------!
  !  Note: The symbols for MU and ALPHA are REVERSED from that of CP2000a,b
  !        Explicit appearance of MUr = 1. has been removed.

  ! Fallspeed parameters:
  real, parameter :: afr=  149.100,  bfr= 0.5000   !Tripoloi and Cotton (1980)
  real, parameter :: afi=   71.340,  bfi= 0.6635   !Ferrier (1994)
  real, parameter :: afs=   11.720,  bfs= 0.4100   !Locatelli and Hobbs (1974)
  real, parameter :: afg=   19.300,  bfg= 0.3700   !Ferrier (1994)
  real, parameter :: afh=  206.890,  bfh= 0.6384   !Ferrier (1994)
 !options:
 !real, parameter :: afs=    8.996,  bfs= 0.4200   !Ferrier (1994)
 !real, parameter :: afg=   6.4800,  bfg= 0.2400   !LH74 (grpl-like snow of lump type)

  real, parameter :: epsQ  = 1.e-14   !min. allowable mixing ratio
  real, parameter :: epsN  = 1.e-3    !min. allowable number concentration
 !real, parameter :: epsZ  = 1.e-32   !min. allowable reflectivity
  real, parameter :: iLAMmin1= 1.e-6  !min. iLAMx (prevents underflow in Nox and VENTx calcs)
  real, parameter :: iLAMmin2= 1.e-10 !min. iLAMx (prevents underflow in Nox and VENTx calcs)
  real, parameter :: eps   = 1.e-32
  real, parameter :: k1    = 0.001
  real, parameter :: k2    = 0.0005
  real, parameter :: k3    = 2.54
  real, parameter :: CPW   = 4218., CPI=2093.
  real, parameter :: deg   =  400., mgo= 1.6e-10
  real, parameter :: deh   =  900.
  real, parameter :: dei   =  500., mio=1.e-12, Nti0=1.e3
  real, parameter :: dew   = 1000.
  real, parameter :: des   =  100., mso=4.4e-10,  rso= 1.0e-4
  real, parameter :: dmr   =    3., dmi=3.,     dms=3.,     dmg=3.,    dmh=3.

  ! NOTE: VxMAX below are the max.allowable mass-weighted fallspeeds for sedimentation.
  !       Thus, Vx corresponds to DxMAX (at sea-level) times the max. density factor, GAM.
  !       [GAMmax=sqrt(DEo/DEmin)=sqrt(1.25/0.4)~2.]  e.g. VrMAX = 2.*8.m/s = 16.m/s
  real, parameter :: DrMAX=  5.e-3,   VrMAX= 16.,   epsQr_sedi= 1.e-8
  real, parameter :: DiMAX= 10.e-3,   ViMAX=  4.,   epsQi_sedi= 1.e-10
  real, parameter :: DsMAX= 30.e-3,   VsMAX=  6.,   epsQs_sedi= 1.e-8
  real, parameter :: DgMAX= 50.e-3,   VgMAX=  8.,   epsQg_sedi= 1.e-8
  real, parameter :: DhMAX= 80.e-3,   VhMAX= 25.,   epsQh_sedi= 1.e-10
  real, parameter :: thrd    = 1./3.
  real, parameter :: sixth   = 0.5*thrd
  real, parameter :: Eci= 1., Ecs= 1., Ecg= 1.
  real, parameter :: Eri= 1., Ers= 1., Erg= 1., Erh= 1.
  real, parameter :: Xdisp   = 0.25               !dispersion of the fall velocity of ice
  real, parameter :: aa11    = 9.44e15, aa22= 5.78e3, Rh= 41.e-6
  real, parameter :: Avx     = 0.78, Bvx= 0.30    !ventilation coefficients [F94 (B.36)]
  real, parameter :: Abigg   = 0.66, Bbigg= 100.  !parameters in probabilistic freezing
  real, parameter :: fdielec = 4.464              !ratio of dielectric factor, |K|w**2/|K|i**2
  real, parameter :: zfact   = 1.e+18             !conversion factor for m-3 to mm2 m-6 for Ze
  real, parameter :: minZET  = -99.               !min. threshold for ZET [dBZ]
  real, parameter :: maxVIS  =  99.e+3            !max. allowable VIS (visibility) [km])
  real, parameter :: Drshed  = 0.001              !mean diam. of drop shed during wet growth
  real, parameter :: SIGcTHRS= 15.e-6             !threshold cld std.dev. before autoconversion
  real, parameter :: KK1= 3.03e3, KK2= 2.59e15    !parameters in Long (1974) kernel
  real, parameter :: Dhh         = 82.e-6         !diameter that rain hump first appears
  real, parameter :: gzMax_sedi  = 180000.        !GZ value below which sedimentation is computed
  real, parameter :: Dr_large    = 200.e-6        !size threshold [m] to distinguish rain/drizzle for precip rates
  real, parameter :: Ds_large    = 200.e-6        !size threshold [m] to distinguish snow/snow-grains for precip rates
  real, parameter :: Dh_large    = 1.0e-2         !size threshold [m] for "large" hail precipitation rate
  real, parameter :: Dr_3cmpThrs = 1.0e-3         !size threshold [m] for hail production from 3-comp freezing
  real, parameter :: Dg_CNgh     = 2.5e-3         !size threshold [m] for CNgh
  real, parameter :: w_CNgh      = 3.             !vertical motion  threshold [m s-1] for CNgh
  real, parameter :: r_CNgh      = 0.05           !Dg/Dho ratio threshold for CNgh
  real, parameter :: Qr_FZrh     = 1.0e-4         !qr-threshold [kg kg-1] for FZrh
  real, parameter :: Tc_FZrh     = -10.           !temp-threshold (C) for FZrh
  real, parameter :: CNsgThres   = 3.0            !threshold for CLcs/VDvs ratio for CNsg
  real, parameter :: capFact     = 0.25           !capacitace of snow
  real, parameter :: qReducFact  = 0.0            !reduction factor for supersaturation water vapor
  real, parameter :: noVal_h_XX  = -1.            !non-value indicator for h_CB, h_ML1, h_ML2, h_SN
  real, parameter :: minSnowSize = 1.e-4          ![m] snow size threshold to compute h_SN

#include "consphy.cdk"
#include "dintern.cdk"
#include "fintern.cdk"

  !Constants used for contact ice nucleation:
  real, parameter :: LAMa0  = 6.6e-8     !mean free path at T0 and p0 [W95_eqn58]
  real, parameter :: T0     = 293.15     !ref. temp.
  real, parameter :: p0     = 101325.    !ref. pres.
  real, parameter :: Ra     = 1.e-6      !aerosol (IN) radius         [M92 p.713; W95_eqn60]
  real, parameter :: kBoltz = 1.381e-23  !Boltzmann's constant
  real, parameter :: KAPa   = 5.39e5     !aerosol thermal conductivity

 !Test switches:
  logical, parameter :: icephase_ON   = .true.  !.false. to suppress ice-phase (Part I)
  logical, parameter :: iceDep_ON     = .true.  !.false. to suppress depositional growth of ice
  logical, parameter :: snow_ON       = .true.  !.false. to suppress snow initiation
  logical, parameter :: grpl_ON       = .true.  !.false. to suppress graupel initiation
  logical, parameter :: hail_ON       = .true.  !.false. to suppress hail initiation
  logical, parameter :: warmphase_ON  = .true.  !.false. to suppress warm-phase (Part II)
  logical, parameter :: autoconv_ON   = .true.  ! autoconversion ON/OFF
  logical, parameter :: rainAccr_ON   = .true.  ! rain accretion and self-collection ON/OFF
  logical, parameter :: sedi_ON       = .true.  !.false. to suppress sedimentation
  logical, parameter :: precipDiag_ON = .true.  !.false. to suppress calc. of sfc precip types
  integer, parameter :: outfreq       =    60   !frequency to compute output diagnostics [s]

  real, dimension(size(QC,dim=1),size(QC,dim=2)) :: DE,DP,QSS,QSW,QSI,WZ,DZ,RHOQX,FLIM,VQQ,  &
        gamfact,gamfact_r,massFlux3D_r,massFlux3D_i,massFlux3D_s,massFlux3D_g
  real, dimension(size(QC,dim=1))                :: HPS,dum
  integer, dimension(size(QC,dim=1))             :: activeColumn

  !==================================================================================!


  !----------------------------------------------------------------------------------!
  !                      PART 1:   Prelimiary Calculations                           !
  !----------------------------------------------------------------------------------!

 !Determine number of upper-levels to exclude from sedimentation:
  ktop_sedi= 0
  do k=1,nk
    ktop_sedi= k
    if (GZ(1,k)<gzMax_sedi) exit
  enddo

 !Compute diagnostic values only every 1 h:
  calcDiag= (mod(nint(dt)*KOUNT,outfreq)==0)

!####  These need only to be computed once per model integration:
!#     (Note:  The values must be available from on all processors)
!#
! if (KOUNT==0) then

   PI2    = PI*2.
   PIov4  = 0.25*PI
   PIov6  = PI*sixth
   CHLS   = CHLC+CHLF  !J k-1; latent heat of sublimation
   LCP    = CHLC/CPD
   LFP    = CHLF/CPD
   iCHLF  = 1./CHLF
   LSP    = LCP+LFP
   ck5    = 4098.170*LCP
   ck6    = 5806.485*LSP
   idew   = 1./dew
   idt    = 1./dt

   !-------------------------------------------------------------------!
   ! Constants based on size distribution parameters:

   ! Mass parameters [ m(D) = cD^d ]
   cmr    = PIov6*dew;  icmr= 1./cmr
   cmi    = 440.;       icmi= 1./cmi
   cms    = PIov6*des;  icms= 1./cms
   cmg    = PIov6*deg;  icmg= 1./cmg
   cmh    = PIov6*deh;  icmh= 1./cmh

   ! Cloud:
   iMUc   =  1./MUc
   GC1    =  gamma(alpha_c+1.0);            iGC1= 1./GC1
   GC2    =  gamma(alpha_c+1.+3.0*iMUc)  !i.e. gamma(alf + 4)
   GC3    =  gamma(alpha_c+1.+6.0*iMUc)  !i.e. gamma(alf + 7)
   GC4    =  gamma(alpha_c+1.+9.0*iMUc)  !i.e. gamma(alf + 10)
   GC11   =  gamma(1.0*iMUc+1.0+alpha_c)
   GC12   =  gamma(2.0*iMUc+1.0+alpha_c)
   GC5    =  gamma(1.0+alpha_c);            iGC5= 1./GC5
   GC6    =  gamma(1.0+alpha_c+1.0*iMUc)
   GC7    =  gamma(1.0+alpha_c+2.0*iMUc)
   GC8    =  gamma(1.0+alpha_c+3.0*iMUc)
   GC13   =  gamma(3.0*iMUc+1.0+alpha_c)
   GC14   =  gamma(4.0*iMUc+1.0+alpha_c)
   GC15   =  gamma(5.0*iMUc+1.0+alpha_c)
   icexc9 =  1./(GC2*iGC1*PIov6*dew)
  !specify cloud droplet number concentration [m-3] based on 'CCNtype':
   if     (CCNtype==1) then
      N_c =  0.8e+8          !maritime
   elseif (CCNtype==2) then
      N_c =  2.0e+8          !continental 1
   elseif (CCNtype==3) then
      N_c =  5.0e+8          !continental 2 (polluted)
   else
      N_c =  2.0e+8          !default (cont1), if 'CCNtype' specified incorrectly
   endif

   ! Rain:
   cexr1  = 1.+alpha_r+dmr+bfr
   cexr2  = 1.+alpha_r+dmr
   ckQr1  = afr*gamma(1.+alpha_r+dmr+bfr)/gamma(1.+alpha_r+dmr)
   GR17   = gamma(2.5+alpha_r+0.5*bfr)
   GR31   =   1. !gamma(1.+alpha_r)
   iGR31  =   1./GR31
   GR32   =   1. !gamma(2.+alpha_r)
   GR33   =   2. !gamma(3.+alpha_r)
   GR34   =   6. !gamma(4.+alpha_r)
   GR35   =  24. !gamma(5.+alpha_r)
   GR36   = 120. !gamma(6.+alpha_r)
   GR37   = 720. !gamma(7.+alpha_r)
   GR50   = (No_r*GR31)**0.75
   cexr5  = 2.+alpha_r
   cexr6  = 2.5+alpha_r+0.5*bfr
   cexr9  = cmr*GR34*iGR31;    icexr9= 1./cexr9
   cexr3  = 1.+bfr+alpha_r
   cexr4  = 1.+alpha_r
   ckQr1  = afr*gamma(1.+dmr+alpha_r+bfr)/gamma(1.+dmr+alpha_r)
   ckQr2  = afr*gamma(1.+bfr+alpha_r)*GR31
   ckQr3  = afr*gamma(7.+alpha_r+bfr)/GR37

   ! Ice:
   GI4    = gamma(alpha_i+dmi+bfi)
   GI6    = gamma(2.5+bfi*0.5+alpha_i)
   GI11   = gamma(1.+bfi+alpha_i)
   GI20   = gamma(0.+bfi+1.+alpha_i)
   GI21   = gamma(1.+bfi+1.+alpha_i)
   GI22   = gamma(2.+bfi+1.+alpha_i)
   GI31   =   1.  !gamma(1.+alpha_i)
   iGI31  =   1./GI31
   GI32   =   1.  !gamma(2.+alpha_i)
   GI33   =   2.  !gamma(3.+alpha_i)
   GI34   =   6.  !gamma(4.+alpha_i)
   GI35   =  24.  !gamma(5.+alpha_i)
   GI36   = 120.  !gamma(6.+alpha_i)
   GI40   = gamma(1.+alpha_i+dmi)
   icexi9 = 1./(cmi*gamma(1.+alpha_i+dmi)*iGI31)
   ckQi1  = afi*gamma(1.+alpha_i+dmi+bfi)/GI40
   ckQi2  = afi*GI11*iGI31
   ckQi4  = 1./(cmi*GI40*iGI31)

   ! Snow:
   cexs1  = 2.5+0.5*bfs+alpha_s
   GS09   = gamma(2.5+bfs*0.5+alpha_s)
   GS11   = gamma(1.+bfs+alpha_s)
   GS12   = gamma(2.+bfs+alpha_s)
   GS13   = gamma(3.+bfs+alpha_s)
   GS31   =   1.  !gamma(1.+alpha_s)
   iGS31  =   1./GS31
   GS32   =   1.  !gamma(2.+alpha_s)
   GS33   =   2.  !gamma(3.+alpha_s)
   GS34   =   6.  !gamma(4.+alpha_s)
   GS35   =  24.  !gamma(5.+alpha_s)
   GS36   = 120.  !gamma(6.+alpha_s)
   GS40   = gamma(1.+alpha_s+dms)
   iGS20  = 1./(GS40*iGS31*cms)
   ckQs1  = afs*gamma(1.+alpha_s+dms+bfs)/GS40
   ckQs2  = afs*GS11*iGS31
   ckQs4  = 1./(cms*GS40*iGS31)

   ! Graupel:
   GG09   = gamma(2.5+0.5*bfg+alpha_g)
   GG11   = gamma(1.+bfg+alpha_g)
   GG12   = gamma(2.+bfg+alpha_g)
   GG13   = gamma(3.+bfg+alpha_g)
   GG31   =   1.  !gamma(1.+alpha_g)
   iGG31  =   1./GG31
   GG32   =   1.  !gamma(2.+alpha_g)
   GG33   =   2.  !gamma(3.+alpha_g)
   GG34   =   6.  !gamma(4.+alpha_g)
   GG35   =  24.  !gamma(5.+alpha_g)
   GG36   = 120.  !gamma(6.+alpha_g)
   GG40   = gamma(1.+alpha_g+dmg)
   iGG99  = 1./(GG40*iGG31*cmg)
   GG50   = (No_g*GG31)**0.75
   ckQg1  = afg*gamma(1.+alpha_g+dmg+bfg)/GG40
   ckQg2  = afg*GG11*iGG31
   ckQg4  = 1./(cmg*GG40*iGG31)

   ! Hail:
   GH09   = gamma(2.5+bfh*0.5+alpha_h)
   GH11   = gamma(1.+bfh+alpha_h)
   GH12   = gamma(2.+bfh+alpha_h)
   GH13   = gamma(3.+bfh+alpha_h)
   GH31   =   1.  !gamma(1.+alpha_h)
   iGH31  =   1./GH31
   GH32   =   1.  !gamma(2.+alpha_h)
   GH33   =   2.  !gamma(3.+alpha_h)
   iGH34  =   1./6. !1./gamma(4.+alpha_h)
   GH40   = gamma(1.+alpha_h+dmh)
   iGH99  = 1./(GH40*iGH31*cmh)
   GH50   = (No_h*GH31)**0.75
   ckQh1  = afh*gamma(1.+alpha_h+dmh+bfh)/GH40
   ckQh2  = afh*GH11*iGH31
   ckQh4  = 1./(cmh*GH40*iGH31)

! endif  !if (KOUNT=0)
!#
!####

!=======================================================================================!

! Compute variables for BLG sedimentation:

  !Compute thickness of layers for sedimentation calcuation:
  tmp1= 1./GRAV
  do k=2,nk
     DZ(:,k)= (GZ(:,k-1)-GZ(:,k))*tmp1
  enddo
  DZ(:,1)= DZ(:,2)

!Max. Courant number for BLG sedimentation (used to determine npassx and dtx):
  CoMAXr = 3.;  CoMAXi = 5.;  CoMAXs = 5.;  CoMAXg = 5.;  CoMAXh = 3.

  !Note: dzMIN is the grid-spacing between the lowest 2 model levels.  The
  !      amount of time-splitting for sedimentation (i.e. the number of times
  !      BLG4 is called each model time step [npassx]) is estimated from this
  !      min DZ, the maximum possible fall speed (VxMAX), and the specified
  !      max Courant number (which may be greater than 1.0 for BLG scheme).

  npassr= max(1, nint( dt*Vrmax/(CoMAXr*dzMIN) ))    !should put in my_mod_sedi, with dzMIN computed
  npassi= max(1, nint( dt*Vimax/(CoMAXi*dzMIN) ))    !for each column (and column-Vxmax estimated)
  npasss= max(1, nint( dt*Vsmax/(CoMAXs*dzMIN) ))
  npassg= max(1, nint( dt*Vgmax/(CoMAXg*dzMIN) ))
  npassh= max(1, nint( dt*Vhmax/(CoMAXh*dzMIN) ))

  dtr = dt/float(npassr);   cr6= GRAV*dtr
  dti = dt/float(npassi);   ci6= GRAV*dti
  dts = dt/float(npasss);   cs6= GRAV*dts
  dtg = dt/float(npassg);   cg6= GRAV*dtg
  dth = dt/float(npassh);   ch6= GRAV*dth
  ck7 = 1./(dt*GRAV)

!=======================================================================================!

! Temporarily store arrays at time (t*) in order to compute (at the end of subroutine)
! the final VXTEND as VXTEND = ( VX{t+1} - VX{t*} )/dt :
  T_TEND = T ;  Q_TEND = Q
  QCTEND = QC;  QRTEND = QR;  QITEND = QI;  QNTEND = QN;  QGTEND = QG;  QHTEND = QH

! Clip all moments to zero if one or more corresponding category moments
!  are less than the minimum allowable value:
!  note: Clipped mass is added back to water vapor field to conserve total mass

  do k= 1,nk
     do i= 1,ni

        if(QC(i,k)<epsQ)    then
           Q(i,k) = Q(i,k) + QC(i,k)
           QC(i,k)= 0.
        endif
        if (QR(i,k)<epsQ)   then
           Q(i,k) = Q(i,k) + QR(i,k)
           QR(i,k)= 0.
        endif
        if (QI(i,k)<epsQ)   then
           Q(i,k) = Q(i,k) + QI(i,k)
           QI(i,k)= 0.
        endif
        if (QN(i,k)<epsQ)   then
           Q(i,k) = Q(i,k) + QN(i,k)
           QN(i,k)= 0.
        endif
        if (QG(i,k)<epsQ)   then
           Q(i,k) = Q(i,k) + QG(i,k)
           QG(i,k)= 0.
        endif
        if (QH(i,k)<epsQ)   then
           Q(i,k) = Q(i,k) + QH(i,k)
           QH(i,k)= 0.
        endif

        if(QCM(i,k)<epsQ)  then
           QM(i,k) = QM(i,k) + QCM(i,k)
           QCM(i,k)= 0.
        endif
        if (QRM(i,k)<epsQ) then
           QM(i,k) = QM(i,k) + QRM(i,k)
           QRM(i,k)= 0.
        endif
        if (QIM(i,k)<epsQ) then
           QM(i,k) = QM(i,k) + QIM(i,k)
           QIM(i,k)= 0.
        endif
        if (QNM(i,k)<epsQ) then
           QM(i,k) = QM(i,k) + QNM(i,k)
           QNM(i,k)= 0.
        endif
        if (QGM(i,k)<epsQ) then
           QM(i,k) = QM(i,k) + QGM(i,k)
           QGM(i,k)= 0.
        endif
        if (QHM(i,k)<epsQ) then
           QM(i,k) = QM(i,k) + QHM(i,k)
           QHM(i,k)= 0.
        endif

    enddo  !i-loop
  enddo    !k-loop;    (clipping)

  QM = max(QM,0.)
  Q  = max(Q ,0.)

! Approximate values at time {t}:
!  [ ave. of values at {*} (advected, but no physics tendency added) and {t-dt} ]:
  HPS= 0.5*(PSM+PS);   TM = 0.5*(TM + T);   QM = 0.5*(QM + Q)
  QCM= 0.5*(QCM+QC);   QRM= 0.5*(QRM+QR);   QIM= 0.5*(QIM+QI)
  QNM= 0.5*(QNM+QN);   QGM= 0.5*(QGM+QG);   QHM= 0.5*(QHM+QH)

  do k=1,nk
     do i=1,ni
    ! Saturation mixing ratios:  [used in both Parts I and II]
        QSW(i,k)= sngl(FOQSA(TM(i,k),HPS(i)*S(i,k)))      !wrt. liquid water at (t)
        QSS(i,k)= sngl(FOQST( T(i,k), PS(i)*S(i,k)))      !wrt. ice surface  at (*)
        QSI(i,k)= sngl(FOQST(TM(i,k),HPS(i)*S(i,k)))      !wrt. ice surface  at (t)
    ! Air density at time (t)
        DE(i,k)= S(i,k)*HPS(i)/(RGASD*TM(i,k))      !air density at time  (t)
     enddo
  enddo

  do i= 1,ni
! Air-density factor: (for fall velocity computations)
     DEo           = DE(i,nk)
     gamfact(i,:)  = sqrt(DEo/(DE(i,:)))
     gamfact_r(i,:)= sqrt( 1./(DE(i,:)))
! Convert 'W_omega' (on thermodynamic levels) to 'w' (on momentum):
     do k= 2,nk-1
        WZ(i,k)= -0.5/(DE(i,k)*GRAV)*(W_omega(i,k-1)+W_omega(i,k+1))
     enddo
     WZ(i,1) = -0.5/(DE(i,1) *GRAV)*W_omega(i,1)
     WZ(i,nk)= -0.5/(DE(i,nk)*GRAV)*W_omega(i,nk)
  enddo

  !----------------------------------------------------------------------------------!
  !                 End of Preliminary Calculation section (Part 1)                  !
  !----------------------------------------------------------------------------------!

  !----------------------------------------------------------------------------------!
  !                      PART 2: Cold Microphysics Processes                         !
  !----------------------------------------------------------------------------------!

! Determine the active grid points (i.e. those which scheme should treat):
  activePoint = .false.
  DO k=2,nk
     DO i=1,ni
        log1= ((QIM(i,k)+QGM(i,k)+QNM(i,k)+QHM(i,k))<epsQ)  !no solid  (i,g,s,h)
        log2= ((QCM(i,k)+QRM(i,k))                  <epsQ)  !no liquid (c,r)
        log3= ((TM(i,k)>TRPL) .and. log1)                   !T>0C & no i,g,s,h
        log4= log1.and.log2.and.(QM(i,k)<QSI(i,k))          !no sol. or liq.; subsat(i)
        if (.not.( log3 .or. log4 ) .and. icephase_ON) then
          activePoint(i,k)= .true.
        endif
     ENDDO
  ENDDO

    ! Size distribution parameters:
    !  Note: + 'thrd' should actually be '1/dmx'(but dmx=3 for all categories x)
    !        + If Qx=0, LAMx etc. are never be used in any calculations
    !          (If Qc=0, CLcy etc. will never be calculated. iLAMx is set to 0
    !           to avoid possible problems due to bugs.)

  DO k= 2,nk                       !Main loop for Part 2
    DO i= 1,ni
      IF (activePoint(i,k)) THEN

    ! Cloud:
       if (QCM(i,k)>epsQ) then
          iQCM   = 1./QCM(i,k)
          iNCM   = 1./N_c
          Dc     = (DE(i,k)*QCM(i,k)*iNCM*icmr)**thrd
          iLAMc  = (DE(i,k)*QCM(i,k)*iNCM*icexc9)**thrd
          iLAMc2 = iLAMc *iLAMc
          iLAMc3 = iLAMc2*iLAMc
          iLAMc4 = iLAMc2*iLAMc2
          iLAMc5 = iLAMc3*iLAMc2
       else
          Dc     = 0.;   iLAMc3= 0.
          iLAMc  = 0.;   iLAMc4= 0.
          iLAMc2 = 0.;   iLAMc5= 0.
       endif

    ! Rain:
       if (QRM(i,k)>epsQ) then
          N_r= GR50*sqrt(sqrt(GR31/GR34*DE(i,k)*QRM(i,k)*icmr))
          iQRM   = 1./QRM(i,k)
          iNRM   = 1./N_r
          Dr     = (DE(i,k)*QRM(i,k)*iNRM*icmr)**thrd
          iLAMr  = max( (DE(i,k)*QRM(i,k)*iNRM*icexr9)**thrd, iLAMmin1 )
          tmp1   = 1./iLAMr
          iLAMr2 = iLAMr *iLAMr
          iLAMr3 = iLAMr2*iLAMr
          iLAMr4 = iLAMr2*iLAMr2
          iLAMr5 = iLAMr3*iLAMr2
          if (Dr>40.e-6) then
             vr0 = gamfact_r(i,k)*ckQr1*iLAMr**bfr
          else
             vr0 = 0.
          endif
       else
          iLAMr  = 0.;  Dr    = 0.;  vr0   = 0.
          iLAMr2 = 0.;  iLAMr3= 0.;  iLAMr4= 0.;  iLAMr5 = 0.
       endif

    ! Ice:
       if (QIM(i,k)>epsQ) then
          N_i    = 5.*exp(0.304*(TRPL-max(233.,TM(i,k))))  !Cooper eqn.
          iQIM   = 1./QIM(i,k)
          iNYM   = 1./N_i
          iLAMi  = max( iLAMmin2, (DE(i,k)*QIM(i,k)*iNYM*icexi9)**thrd )
          iLAMi2 = iLAMi *iLAMi
          iLAMi3 = iLAMi2*iLAMi
          iLAMi4 = iLAMi2*iLAMi2
          iLAMi5 = iLAMi3*iLAMi2
          iLAMiB0= iLAMi**(bfi)
          iLAMiB1= iLAMi**(bfi+1.)
          iLAMiB2= iLAMi**(bfi+2.)
          vi0    = gamfact(i,k)*ckQi1*iLAMiB0
       else
          iLAMi  = 0.;  vi0    = 0.
          iLAMi2 = 0.;  iLAMi3 = 0.;  iLAMi4 = 0.;  iLAMi5= 0.
          iLAMiB0= 0.;  iLAMiB1= 0.;  iLAMiB2= 0.
       endif

    ! Snow:
       if (QNM(i,k)>epsQ) then
          No_s   = min(2.e+8, 2.e+6*exp(-0.12*min(-0.001,TM(i,k)-TRPL))) !T2004
          N_s    = (No_s*GS31)**0.75*sqrt(sqrt(GS31/GS34*DE(i,k)*QNM(i,k)*icms))
          iQNM   = 1./QNM(i,k)
          iNNM   = 1./N_s
          iLAMs  = max(iLAMmin2, (DE(i,k)*QNM(i,k)*iNNM*iGS20)**thrd )
          iLAMs2 = iLAMs*iLAMs
          iLAMsB0= iLAMs**(bfs)
          iLAMsB1= iLAMs**(bfs+1.)
          iLAMsB2= iLAMs**(bfs+2.)
          vs0    = gamfact(i,k)*ckQs1*iLAMsB0
       else
          iLAMs  = 0.;  vs0    = 0.;  iLAMs2= 0.
          iLAMsB0= 0.;  iLAMsB1= 0.;  iLAMsB1= 0.
       endif

    ! Graupel:
       if (QGM(i,k)>epsQ) then
          N_g    = GG50*sqrt(sqrt(GG31/GG34*DE(i,k)*QGM(i,k)*icmg))
          iQGM   = 1./QGM(i,k)
          iNGM   = 1./N_g
          iLAMg  = max(iLAMmin1, (DE(i,k)*QGM(i,k)*iNGM*iGG99)**thrd )
          iLAMg2 = iLAMg *iLAMg
          iLAMgB0= iLAMg**(bfg)
          iLAMgB1= iLAMg**(bfg+1.)
          iLAMgB2= iLAMg**(bfg+2.)
          vg0    = gamfact(i,k)*ckQg1*iLAMgB0
       else
          iLAMg  = 0.;  vg0  = 0.
          iLAMg2 = 0.;  iLAMgB0= 0.;  iLAMgB1= 0.;  iLAMgB1= 0.
       endif

    ! Hail:
       if (QHM(i,k)>epsQ) then
          N_h    = GH50*sqrt(sqrt(GH31*iGH34*DE(i,k)*QHM(i,k)*icmh))
          iQHM   = 1./QHM(i,k)
          iNHM   = 1./N_h
          iLAMh  = max(iLAMmin1, (DE(i,k)*QHM(i,k)*iNHM*iGH99)**thrd)
          iLAMh2 = iLAMh*iLAMh
          iLAMhB0= iLAMh**(bfh)
          iLAMhB1= iLAMh**(bfh+1.)
          iLAMhB2= iLAMh**(bfh+2.)
          vh0    = gamfact(i,k)*ckQh1*iLAMhB0
       else
          iLAMh  = 0.;  vh0    = 0.
          iLAMhB0= 0.;  iLAMhB1= 0.;  iLAMhB1= 0.
       endif
!------

 !Calculating ice-phase source/sink terms:

 ! Initialize all source terms to zero:
       QNUvi=0.;  QVDvi=0.;  QVDvs=0.;  QVDvg=0.;  QVDvh=0.;   QCLci=0.
       QCLcs=0.;  QCLcg=0.;  QCLch=0.;  QFZci=0.;  QCLri=0.;   QMLsr=0.
       QCLrs=0.;  QCLrg=0.;  QMLgr=0.;  QCLrh=0.;  QMLhr=0.;   QFZrh=0.
       QMLir=0.;  QCLci=0.;  QCNig=0.;  QCLsr=0.;  QCLsh=0.;   QCLgr=0.
       QCNis=0.;  QCLir=0.;  QCLis=0.;  QCLig=0.;  QCLih=0.;   QCNgh=0.
       QIMsi=0.;  QIMgi=0.;  QCNsg=0.;  QHwet=0.
       Dirg =0.;  Dirh=0.; Dsrs= 0.; Dsrg= 0.; Dsrh= 0.; Dgrg=0.; Dgrh=0.

       Tc    = TM(i,k)-TRPL
       if (Tc<-120.)   &
        print*, '***WARNING*** -- In MULTIMOMENT --  Ambient Temp.(C):',Tc
       Cdiff = (2.2157e-5+0.0155e-5*Tc)*1.e5/(S(i,k)*HPS(i))
       MUdyn = 1.72e-5*(393./((TM(i,k)+120.)))*(TM(i,k)/TRPL)**1.5 !RYp.102
       MUkin = MUdyn/(DE(i,k))
       iMUkin= 1./MUkin
       ScTHRD= (MUkin/Cdiff)**thrd       ! i.e. Sc^(1/3)
       Ka    = 2.3971e-2 + 0.0078e-2*Tc                                   !therm.cond.(air)
       Kdiff = (9.1018e-11*TM(i,k)*TM(i,k)+8.8197e-8*TM(i,k)-(1.0654e-5)) !therm.diff.(air)
       iDE   = 1./DE(i,k)
       gam   = gamfact(i,k)

        !Collection efficiencies:
       Eis   = min(0.05*exp(0.1*Tc),1.)     !F95 (Table 1)
       Eig   = min(0.01*exp(0.1*Tc),1.)     !dry (Eig=1.0 for wet growth)
       Eii   = 0.1*Eis
       Ess   = Eis;   Eih = Eig;   Esh = Eig
       iEih  = 1./Eih
       iEsh  = 1./Esh
       !NOTE:  Eci=Ecs=Ecg=Eri=Ers=Erg=Erh=1. (constant parameters)
       !       Ech is computed in CLch section

       qvs0  = sngl(FOQSA(TRPL,HPS(i)*S(i,k)))      !sat.mix.ratio at 0C
       DELqvs= qvs0-(QM(i,k))

   !-------------------------------------------------------------------------------------------!

           ! COLLECTION by snow, graupel, hail:
           !  (i.e. wet or dry ice-categories [=> excludes ice crystals])

           ! Collection by SNOW:
       if (QNM(i,k)>epsQ) then
          ! cloud:
          if (QCM(i,k)>epsQ) then
             QCLcs= dt*gam*afs*cmr*Ecs*PIov4*iDE*(N_c*N_s)*iGC5*iGS31*(GC13*GS13*       &
                    iLAMc3*iLAMsB2+2.*GC14*GS12*iLAMc4*iLAMsB1+GC15*GS11*iLAMc5*        &
                    iLAMsB0)
             QCLcs= min(QCLcs, (QCM(i,k)))
          else
             QCLcs= 0.;   NCLcs= 0.
          endif

          ! ice:
          if (QIM(i,k)>epsQ) then
             tmp1 = vs0-vi0
             tmp3 = sqrt(tmp1*tmp1+0.04*vs0*vi0)
             QCLis= dt*cmi*iDE*PI*6.*Eis*(N_i*N_s)*tmp3*iGI31*iGS31*(0.5*iLAMs2*iLAMi3+ &
                    2.*iLAMs*iLAMi4+5.*iLAMi5)
             QCLis= min(QCLis, (QIM(i,k)))
          else
             QCLis= 0.
          endif

       else
          QCLcs= 0.;  QCLis= 0.
       endif

       ! Collection by GRAUPEL:
       if (QGM(i,k)>epsQ) then

          ! cloud:
          if (QCM(i,k)>epsQ) then
             QCLcg= dt*gam*afg*cmr*Ecg*PIov4*iDE*(N_c*N_g)*iGC5*iGG31*(GC13*GG13*       &
                    iLAMc3*iLAMgB2+ 2.*GC14*GG12*iLAMc4*iLAMgB1+GC15*GG11*iLAMc5*       &
                    iLAMgB0)
             QCLcg= min(QCLcg, (QCM(i,k)))
          else
             QCLcg= 0.
          endif

          ! ice:
          if (QIM(i,k)>epsQ) then
             tmp1 = vg0-vi0
             tmp3 = sqrt(tmp1*tmp1+0.04*vg0*vi0)
             QCLig= dt*cmi*iDE*PI*6.*Eig*(N_i*N_g)*tmp3*iGI31*iGG31*(0.5*iLAMg2*iLAMi3+ &
                    2.*iLAMg*iLAMi4+5.*iLAMi5)
             QCLig= min(QCLig, (QIM(i,k)))
          else
             QCLig= 0.
          endif

       else
          QCLcg= 0.;   QCLrg= 0.;   QCLig= 0.
       endif

       ! Collection by HAIL:
       if (QHM(i,k)>epsQ) then

         ! cloud:
          if (QCM(i,k)>epsQ) then
             Dh   = ((DE(i,k)*QHM(i,k)*iNHM)*icmh)**thrd
             Ech  = exp(-8.68e-7*Dc**(-1.6)*Dh)    !Z85_A24
             QCLch= dt*gam*afh*cmr*Ech*PIov4*iDE*(N_c*N_h)*iGC5*iGH31*(GC13*GH13*       &
                    iLAMc3*iLAMhB2+2.*GC14*GH12*iLAMc4*iLAMhB1+GC15*GH11*iLAMc5*iLAMhB0)
             QCLch= min(QCLch, QCM(i,k))
          else
             QCLch= 0.
          endif

          ! rain:
          if (QRM(i,k)>epsQ) then
             tmp1 = vh0-vr0
             tmp3 = sqrt(tmp1*tmp1+0.04*vh0*vr0)
             QCLrh= dt*cmr*Erh*PIov4*iDE*(N_h*N_r)*iGR31*iGH31*tmp3*(GR36*GH31*iLAMr5+  &
                    2.*GR35*GH32*iLAMr4*iLAMh+GR34*GH33*iLAMr3*iLAMh2)
             QCLrh= min(QCLrh, QRM(i,k))
          else
             QCLrh= 0.
          endif

          ! ice:
          if (QIM(i,k)>epsQ) then
             tmp1 = vh0-vi0
             tmp3 = sqrt(tmp1*tmp1+0.04*vh0*vi0)
             QCLih= dt*cmi*iDE*PI*6.*Eih*(N_i*N_h)*tmp3*iGI31*iGH31*(0.5*iLAMh2*iLAMi3+ &
                    2.*iLAMh*iLAMi4+5.*iLAMi5)
             QCLih= min(QCLih, QIM(i,k))
          else
             QCLih= 0.
          endif

          ! snow:
          if (QNM(i,k)>epsQ) then
             tmp1 = vh0-vs0
             tmp3 = sqrt(tmp1*tmp1+0.04*vh0*vs0)
             tmp4 = iLAMs2*iLAMs2
             QCLsh= dt*cms*iDE*PI*6.*Esh*(N_s*N_h)*tmp3*iGS31*iGH31*(0.5*iLAMh2*iLAMs2* &
                    iLAMs+2.*iLAMh*tmp4+5.*tmp4*iLAMs)
             QCLsh= min(QCLsh, (QNM(i,k)))
          else
             QCLsh= 0.
          endif

          ! wet growth:
          VENTh= Avx*GH32*iLAMh*iLAMh + Bvx*ScTHRD*sqrt(gam*afh*iMUkin)*GH09*iLAMh**    &
                 (2.5+0.5*bfh+alpha_h)
          QHwet= max(0., dt*PI2*(DE(i,k)*CHLC*Cdiff*DELqvs-Ka*Tc)*No_h*iDE/(CHLF+CPW*   &
                 Tc)*VENTh+(QCLih*iEih+QCLsh*iEsh)*(1.-CPI*Tc/(CHLF+CPW*Tc)) )

       else
          QCLch= 0.;   QCLrh= 0.;   QCLih= 0.;   QCLsh= 0.;   QHwet= 0.
       endif

       IF (TM(i,k)>TRPL .and. warmphase_ON) THEN
          !**********!
          !  T > To  !
          !**********!

          ! MELTING of frozen particles:
         !  ICE:
          QMLir   = (QIM(i,k))
          QIM(i,k)= 0.
          NMLir   = (N_i)

          !  SNOW:
          if (QNM(i,k)>epsQ) then
             VENTs= Avx*GS32*iLAMs*iLAMs + Bvx*ScTHRD*sqrt(gam*afs*iMUkin)*GS09*iLAMs**  &
                    cexs1
             QMLsr= dt*(PI2*iDE*iCHLF*No_s*VENTs*(Ka*Tc-CHLC*Cdiff*DELqvs) + CPW*iCHLF*  &
                    Tc*(QCLcs+QCLrs)*idt)
             QMLsr= min(max(QMLsr,0.), QNM(i,k))
          else
             QMLsr= 0.
          endif

          !  GRAUPEL:
          if (QGM(i,k)>epsQ) then
             VENTg= Avx*GG32*iLAMg*iLAMg + Bvx*ScTHRD*sqrt(gam*afg*iMUkin)*GG09*iLAMg**  &
                    (2.5+0.5*bfg+alpha_g)
             QMLgr= dt*(PI2*iDE*iCHLF*No_g*VENTg*(Ka*Tc-CHLC*Cdiff*DELqvs) + CPW*iCHLF*  &
                    Tc*(QCLcg+QCLrg)*idt)
             QMLgr= min(max(QMLgr,0.), QGM(i,k))
          else
             QMLgr= 0.
          endif

          !  HAIL:
          if (QHM(i,k)>epsQ.and.Tc>5.) then
             VENTh= Avx*GH32*iLAMh*iLAMh + Bvx*ScTHRD*sqrt(gam*afh*iMUkin)*GH09*iLAMh**  &
                    (2.5+0.5*bfh+alpha_h)
             QMLhr= dt*(PI2*iDE*iCHLF*No_h*VENTh*(Ka*Tc-CHLC*Cdiff*DELqvs) + CPW/CHLF*   &
                    Tc*(QCLch+QCLrh)*idt)
             QMLhr= min(max(QMLhr,0.), QHM(i,k))
          else
             QMLhr= 0.
          endif

         ! Cold (sub-zero) source/sink terms:
          QNUvi= 0.;   QFZci= 0.;   QVDvi= 0.;   QVDvs= 0.;   QVDvg= 0.
          QCLci= 0.;   QCLis= 0.;   QCNig= 0.;   QCNis1=0.;   QCNis2=0.
          QCNgh= 0.;   QIMsi= 0.;   QIMgi= 0.;   QCLir= 0.;   QCLri= 0.
          QCLrs= 0.;   QCLgr= 0.;   QCLrg= 0.;   QCNis= 0.;   QVDvh= 0.
          QCNsg= 0.;   QCLsr= 0.

       ELSE
          !----------!
          !  T < To  !
          !----------!
          tmp1  = 1./QSI(i,k)
          Si    = QM(i,k) *tmp1
          tmp2= TM(i,k)*TM(i,k)
          iABi  = 1./( CHLF*CHLF/(Ka*RGASV*tmp2) + 1./(DE(i,k)*(QSI(i,k))*Cdiff) )

          ! Warm-air-only source/sink terms:
          QMLir= 0.;   QMLsr= 0.;   QMLgr= 0.;   QMLhr= 0.

          ! Probabilistic freezing (Bigg) of rain:
          if (QRM(i,k)>Qr_FZrh .and. Tc<Tc_FZrh .and. hail_ON) then
             NrFZrh= -dt*Bbigg*(exp(Abigg*Tc)-1.)*DE(i,k)*QRM(i,k)*idew
             Rz= 1.  !N and Z (and Q) are conserved for FZrh with triple-moment
           ! The Rz factor serves to conserve reflectivity when a rain distribution
           !  converts to an distribution with a different shape parameter, alpha.
           !  (e.g. when rain freezes to hail)  The factor Rz non-conserves N while
           !  acting to conserve Z for double-moment.  See Ferrier, 1994 App. D)
             NhFZrh= Rz*NrFZrh
             QFZrh = NrFZrh*(QRM(i,k)*iNRM)
          else
             QFZrh= 0.
          endif

          !--------!
          !  ICE:  !
          !--------!
          ! Homogeneous freezing of cloud to ice:  (simplified)
          if (QCM(i,k)>epsQ .and. Tc<-35.) then
             FRAC= 1.  !if T<-35
             QFZci= FRAC*QCM(i,k)
          else
             QFZci= 0.
          endif

!  Ice initiation (single-moment):
!         if (QIM(i,k)<=epsQ .and. Tc<-5. .and. Si>1.) then
          if (QIM(i,k)<=epsQ .and. Tc<-5. .and. Si>1.08) then !following Thompson etal (2006)
             NNUvi = 5.*exp(0.304*(TRPL-max(233.,TM(i,k))))  !Cooper eqn.
             QNUvi= mio*iDE*NNUvi
             QNUvi= min(QNUvi,Q(i,k))
          endif

          IF (QIM(i,k)>epsQ) THEN

             ! Riming (stochastic collection of cloud):
             if (QCM(i,k)>epsQ) then
                QCLci= dt*gam*afi*cmr*Eci*PIov4*iDE*(N_c*N_i)*iGC5*iGI31*(GC13*GI22*     &
                       iLAMc3*iLAMiB2+2.*GC14*GI21*iLAMc4*iLAMiB1+GC15*GI20*iLAMc5*      &
                       iLAMiB0)
                QCLci= min(QCLci, (QCM(i,k)))
             else
                QCLci= 0.
             endif

             ! Deposition/sublimation:
             No_i = N_i*iGI31/iLAMi
             VENTi= Avx*GI32*iLAMi*iLAMi + Bvx*ScTHRD*sqrt(gam*afi*iMUkin)*GI6*iLAMi**   &
                    (2.5+0.5*bfi+alpha_i)
             QVDvi= dt*iABi*(PI2*(Si-1.)*No_i*VENTi - idt*QCLci*CHLS*CHLF/(Ka*RGASV*     &
                    TM(i,k)*TM(i,k)))

             ! Prevent overdepletion of vapor:
             tmp1  = T(i,k)-7.66
             VDmax = (Q(i,k)-QSS(i,k))/(1.+ck6*(QSS(i,k))/(tmp1*tmp1))
             if(Si>=1.) then
                QVDvi= min(max(QVDvi,0.),VDmax)
             else
                if (VDmax<0.) QVDvi= max(QVDvi,VDmax)
               !IF prevents subl.(QVDvi<0 at t) changing to dep.(VDmax>0 at t*)  2005-06-28
             endif

             if (.not. iceDep_ON) QVDvi= 0. !to suppress depositional growth

             ! Conversion to graupel:
             QCNig= max(QCLci-max(0.,QVDvi), 0.)
             if (.not. grpl_ON) QCNig= 0.
             QCLci= max(QCLci-QCNig, 0.)

             ! Conversion to snow:
             !   +depostion/riming of ice:
             Di= (DE(i,k)*QIM(i,k)*iNYM*icmi)**thrd
             mi= DE(i,k)*(QIM(i,k)*iNYM)
             ri= 0.5*Di
             if (mi<=0.5*mso.and.abs(0.5*mso-mi)>1.e-20) then
                QCNis1= (mi/(mso-mi))*(QVDvi+QCLci)
             else
                QCNis1= (QVDvi+QCLci) + (1.-0.5*mso/mi)*QIM(i,k)
             endif
             QCNis1= max(0., QCNis1)
             !   +aggregation of ice:
             if(ri<0.5*rso) then
                Ki    = PIov6*Di*Di*vi0*Eii*Xdisp
                tmp1  = log(ri/rso)
                tmp2  = tmp1*tmp1*tmp1
                QCNis2= -dt*0.5*(QIM(i,k)*N_i)*Ki/tmp2
             else
                Ki= 0.;   QCNis2= 0.
             endif
             !   +total conversion rate:
             QCNis =  QCNis1 + QCNis2
             if (.not.(snow_ON)) then
                QCNis = 0.
             endif

             ! 3-component freezing (collisions with rain):
             if (QRM(i,k)>epsQ .and. QIM(i,k)>epsQ) then
                tmp1= vr0-vi0
                tmp2= tmp1*tmp1
                tmp3= sqrt(tmp2+0.04*vr0*vi0)

                QCLir= dt*cmi*Eri*PIov4*iDE*(N_r*N_i)*iGI31*iGR31*tmp3*(GI36*GR31*      &
                       iLAMi5+2.*GI35*GR32*iLAMi4*iLAMr+GI34*GR33*iLAMi3*iLAMr2)
                QCLri= dt*cmr*Eri*PIov4*iDE*(N_i*N_r)*iGR31*iGI31*tmp3*(GR36*GI31*      &
                       iLAMr5+2.*GR35*GI32*iLAMr4*iLAMi+GR34*GI33*iLAMr3*iLAMi2)
                QCLri= min(QCLri, (QRM(i,k)));  QCLir= min(QCLir, (QIM(i,k)))

                !Determine destination of 3-comp.freezing:
                tmp1= max(Di,Dr)  !OPTIM NOTE:  Make sure Di gets calculated
                dey= (dei*Di*Di*Di+dew*Dr*Dr*Dr)/(tmp1*tmp1*tmp1)
                if (dey>0.5*(deg+deh) .and. Dr>Dr_3cmpThrs .and. hail_ON) then
                   Dirg= 0.;  Dirh= 1.
                else
                   Dirg= 1.;  Dirh= 0.
                endif
             else
                QCLir= 0.;  QCLri= 0.
                Dirh = 0.;  Dirg = 0.
             endif
             if (.not. grpl_ON) Dirg= 0.

             !  Rime-splintering (ice multiplication):
             ff= 0.
             if(Tc>=-8..and.Tc<=-5.) ff= 3.5e8*(Tc +8.)*thrd
             if(Tc> -5..and.Tc< -3.) ff= 3.5e8*(-3.-Tc)*0.5
             NIMii= DE(i,k)*ff*QCLci
             NIMsi= DE(i,k)*ff*QCLcs
             NIMgi= DE(i,k)*ff*QCLcg
             QIMsi= mio*iDE*NIMsi
             QIMgi= mio*iDE*NIMgi

          ELSE

             QCLci= 0.;  QVDvi= 0.;   QCNig= 0.;   QCNis= 0.
             QIMsi= 0.;  QIMgi= 0.;   QCLri= 0.;   QCLir= 0.

          ENDIF
          !---------!
          !  SNOW:  !
          !---------!
          IF (QNM(i,k)>epsQ) THEN

             ! Deposition/sublimation:
             VENTs = Avx*GS32*iLAMs*iLAMs + Bvx*ScTHRD*sqrt(gam*afs*iMUkin)*GS09*iLAMs** &
                     (2.5+0.5*bfs+alpha_s)
             QVDvs = dt*capFact*iABi*(PI2*(Si-1.)*No_s*VENTs - CHLS*CHLF/(Ka*RGASV*      &
                     TM(i,k)*TM(i,k))*QCLcs*idt)
             ! Prevent overdepletion of vapor:
             tmp1  = T(i,k)-7.66
             VDmax = (Q(i,k)-QSS(i,k))/(1.+ck6*(QSS(i,k))/(tmp1*tmp1))  !KY97_A.33

             if(Si>=1.) then
                QVDvs= min(max(QVDvs,0.),VDmax)
             else
          !     QVDvs= max(QVDvs,VDmax)
                if (VDmax<0.) QVDvs= max(QVDvs,VDmax)
                !IF prevents subl.(QVDvs<0 at t) changing to dep.(VDmax>0 at t*)  2005-06-28
             endif
             NVDvs= -min(0.,(N_s*iQNM)*QVDvs)  !pos. quantity

             ! Conversion to graupel:
             if (QCLcs > CNsgThres*QVDvs) then
                QCNsg= (deg/(deg-des))*QCLcs
             else
                QCNsg= 0.
             endif
             if (.not. grpl_ON) QCNsg= 0.

             ! 3-component freezing (collisions with rain):
              if (QRM(i,k)>epsQ .and. QNM(i,k)>epsQ .and. Tc<-5.) then
                tmp1 = vs0-vr0
                tmp2 = sqrt(tmp1*tmp1+0.04*vs0*vr0)
                tmp6 = iLAMs2*iLAMs2*iLAMs
                QCLrs= dt*cmr*Ers*PIov4*iDE*N_s*N_r*iGR31*iGS31*tmp2*(GR36*GS31*iLAMr5+ &
                       2.*GR35*GS32*iLAMr4*iLAMs+GR34*GS33*iLAMr3*iLAMs2)
                QCLsr= dt*cms*Ers*PIov4*iDE*(N_r*N_s)*iGS31*iGR31*tmp2*(GS36*GR31*tmp6+ &
                       2.*GS35*GR32*iLAMs2*iLAMs2*iLAMr+GS34*GR33*iLAMs2*iLAMs*iLAMr2)
               !note: For explicit eqns, NCLsr = NCLrs
                QCLrs= min(QCLrs, QRM(i,k));  QCLsr= min(QCLsr, QIM(i,k))

                ! Determine destination of 3-comp.freezing:
                Dsrs= 0.;   Dsrg= 0.;    Dsrh= 0.
                Ds  = (DE(i,k)*QNM(i,k)*iNNM*icms)**thrd
                tmp1= max(Ds,Dr)
                tmp2= tmp1*tmp1*tmp1
                dey = (des*Ds*Ds*Ds + dew*Dr*Dr*Dr)/tmp2
                if (dey<=0.5*(des+deg)                        ) Dsrs= 1.  !snow
                if (dey >0.5*(des+deg) .and. dey<0.5*(deg+deh)) Dsrg= 1.  !graupel
                if (dey>=0.5*(deg+deh)) then
                   Dsrh= 1.                                               !hail
                   if (.not.hail_ON .or. Dr<Dr_3cmpThrs) then
                      Dsrg= 1.;   Dsrh= 0.                                !graupel
                   endif
                endif
                if (.not. grpl_ON) Dsrg= 0.

             else
                QCLrs= 0.;   QCLsr= 0.
             endif

          ELSE

             QVDvs= 0.;  QCLcs= 0.;  QCNsg= 0.;  QCLsr= 0.;  QCLrs= 0.

          ENDIF
          !------------!
          !  GRAUPEL:  !
          !------------!
          IF (QGM(i,k)>epsQ) THEN
             Dg   = ((DE(i,k)*QGM(i,k)*iNGM)*icmg)**thrd

           !Conversion to hail:    (Dho given by S-L limit)
             if (Dg>Dg_CNgh .and. WZ(i,k)>w_CNgh .and. hail_ON) then
          !  note:  The Dg threshold for CNgh to occur is a surrogate for the
          !         presense of large graupel (that should convert to hail under riming).
          !         However, while this is meaningful for single-moment, it is not appropriate
          !         for double-moment, since a given Dg may have a large or small quantity
          !         of "large" graupel.  For d-m, need to change condition based on estimate
          !         of the amount of large graupel (incomplete gamma distribution).
                Dho  = 0.01*(exp(min(20.,-Tc/(1.1e4*DE(i,k)*(QCM(i,k)+QRM(i,k))-1.3e3*   &
                      DE(i,k)*(QIM(i,k))+1.)))-1.)
                Dho  = min(1., max(0.0001,Dho))    !smallest Dho=0.1mm; largest=1m
                ratio= Dg/Dho
                if (ratio>r_CNgh) then
                   QCNgh= (0.5*ratio)*(QCLcg+QCLrg+QCLig)
                   QCNgh= min(QCNgh,(QGM(i,k))+QCLcg+QCLrg+QCLig)
                else
                   QCNgh= 0.
                endif
             endif

             ! 3-component freezing (collisions with rain)
             if (QRM(i,k)>epsQ) then
                tmp1 = vg0-vr0
                tmp2 = sqrt(tmp1*tmp1 + 0.04*vg0*vr0)
                tmp8 = iLAMg2*iLAMg      ! iLAMg**3
                tmp9 = tmp8*iLAMg        ! iLAMg**4
                tmp10= tmp9*iLAMg        ! iLAMg**5
                QCLrg= dt*cmr*Erg*PIov4*iDE*(N_g*N_r)*iGR31*iGG31*tmp2*(GR36*GG31*      &
                       iLAMr5+2.*GR35*GG32*iLAMr4*iLAMg+GR34*GG33*iLAMr3*iLAMg2)
                QCLgr= dt*cms*Erg*PIov4*iDE*(N_r*N_g)*iGG31*iGR31*tmp2*(GG36*GR31*tmp10 &
                       +2.*GG35*GR32*tmp9*iLAMr+GG34*GR33*tmp8*iLAMr2)
               !(note: For explicit eqns, NCLgr= NCLrg)
                QCLrg= min(QCLrg, QRM(i,k));  QCLgr= min(QCLgr, QGM(i,k))

               ! Determine destination of 3-comp.freezing:
                tmp1= max(Dg,Dr)
                tmp2= tmp1*tmp1*tmp1
                dey = (deg*Dg*Dg*Dg + dew*Dr*Dr*Dr)/tmp2
                if (dey>0.5*(deg+deh) .and. Dr>Dr_3cmpThrs .and. hail_ON) then
                   Dgrg= 0.;  Dgrh= 1.
                else
                   Dgrg= 1.;  Dgrh= 0.
                endif
             else
                QCLgr= 0.;  QCLrg= 0.;  NCLgr= 0.;  NCLrg= 0.
             endif

          ELSE

             QVDvg= 0.;  QCNgh= 0.;  QCLgr= 0.;  QCLrg= 0.

          ENDIF
          !---------!
          !  HAIL:  !
          !---------!
          IF (QHM(i,k)>epsQ) THEN
             ! Wet growth:
             if (QHwet<(QCLch+QCLrh+QCLih+QCLsh) .and. Tc>-40.) then
                QCLih= min(QCLih*iEih, QIM(i,k))   !change Eih to 1. in CLih
                QCLsh= min(QCLsh*iEsh, QNM(i,k))   !change Esh to 1. in CLsh
                tmp3 = QCLrh
                QCLrh= QHwet-(QCLch+QCLih+QCLsh)   !actual QCLrh minus QSHhr
                QSHhr= tmp3-QCLrh                  !QSHhr used here only
             endif
          ELSE
             QVDvh= 0.
          ENDIF


       ENDIF  ! ( if Tc<0C Block )

     !------------  End of source/sink term calculation  -------------!

     !Iterative adjustment of sink (and source) terms to prevent overdepletion:
       do niter= 1,2

          ! (1) Vapor:
          source= (Q(i,k)) +dim(-QVDvi,0.)+dim(-QVDvs,0.)+dim(-QVDvg,0.)+dim(-QVDvh,0.)
          sink  = QNUvi+dim(QVDvi,0.)+dim(QVDvs,0.)
          sour  = max(source,0.)
          if(sink>sour) then
             ratio= sour/sink
             QNUvi= ratio*QNUvi
             if(QVDvi>0.) then
               QVDvi= ratio*QVDvi
             endif
             if(QVDvs>0.) then
               QVDvs=ratio*QVDvs
             endif
             QVDvg= ratio*QVDvg;   QVDvh= ratio*QVDvh
          endif

          ! (2) Cloud:
          source= (QC(i,k))
          sink  = QCLci+QCLcs+QCLcg+QCLch+QFZci
          sour  = max(source,0.)
          if(sink>sour) then
             ratio= sour/sink
             QFZci= ratio*QFZci;   QCLci= ratio*QCLci
             QCLcs= ratio*QCLcs;
             QCLcg= ratio*QCLcg;   QCLch= ratio*QCLch;
          endif

          ! (3) Rain:
          source= (QR(i,k))+QMLsr+QMLgr+QMLhr+QMLir
          sink  = QCLri+QCLrs+QCLrg+QCLrh+QFZrh
          sour  = max(source,0.)
          if(sink>sour) then
             ratio= sour/sink
             QCLrg= ratio*QCLrg;   QCLri= ratio*QCLri
             QFZrh= ratio*QFZrh;   QCLrg= ratio*QCLrg
             QCLrs= ratio*QCLrs;   QCLrh= ratio*QCLrh
             if (ratio==0.) then
                Dirg= 0.; Dirh= 0.; Dgrg= 0.; Dgrh= 0.
                Dsrs= 0.; Dsrg= 0.; Dsrh= 0.
              endif
          endif

          ! (4) Ice:
          source= (QI(i,k))+QNUvi+dim(QVDvi,0.)+QCLci+QFZci
          sink  = QCNig+QCNis+QCLir+dim(-QVDvi,0.)+QCLis+QCLig+QCLih+QMLir
          sour  = max(source,0.)
          if(sink>sour) then
             ratio= sour/sink
             QMLir= ratio*QMLir;    NMLir= ratio*NMLir
             if (QVDvi<0.) then
                QVDvi= ratio*QVDvi; NVDvi= ratio*NVDvi
             endif
             QCLig=  ratio*QCLig
             QCNig=  ratio*QCNig
             QCNis=  ratio*QCNis
             QCLir=  ratio*QCLir
             QCLis=  ratio*QCLis
             QCLih=  ratio*QCLih
             if (ratio==0.) then
                Dirg= 0.; Dirh= 0.
             endif
          endif

          ! (5) Snow:
          source= (QN(i,k))+QCNis+dim(QVDvs,0.)+QCLis+Dsrs*(QCLrs+QCLsr)+QCLcs
          sink  = dim(-QVDvs,0.)+QCNsg+QMLsr+QCLsr+QCLsh
          sour  = max(source,0.)
          if(sink>sour) then
             ratio= sour/sink
             if(QVDvs<=0.) then
                QVDvs= ratio*QVDvs;   NVDvs= ratio*NVDvs
             endif
             QCNsg= ratio*QCNsg;   QMLsr= ratio*QMLsr
             QCLsr= ratio*QCLsr;   QCLsh= ratio*QCLsh
             if (ratio==0.) then
                Dsrs= 0.; Dsrg= 0.; Dsrh= 0.
             endif
          endif

          !  (6) Graupel:
          source= (QG(i,k))+QCNig+QCNsg+dim(QVDvg,0.)+Dirg*(QCLri+QCLir)+Dgrg*(QCLrg+    &
                  QCLgr)+QCLcg+Dsrg*(QCLrs+QCLsr)+QCLig
          sink  = dim(-QVDvg,0.)+QMLgr+QCNgh+QCLgr
          sour  = max(source,0.)
          if(sink>sour) then
             ratio= sour/sink
             QVDvg= ratio*QVDvg;   QMLgr= ratio*QMLgr
             QCNgh= ratio*QCNgh;   QCLgr= ratio*QCLgr
             if (ratio==0.) then
                Dgrg= 0.; Dgrh= 0.
             endif
          endif

          !  (7) Hail:
          source= (QH(i,k))+dim(QVDvh,0.)+QCLch+QCLrh+Dirh*(QCLri+QCLir)+QCLih+QCLsh+    &
                  Dsrh*(QCLrs+QCLsr)+QCNgh+Dgrh*(QCLrg+QCLgr)+QFZrh
          sink  = dim(-QVDvh,0.)+QMLhr
          sour  = max(source,0.)
          if(sink>sour) then
             ratio= sour/sink
             QVDvh= ratio*QVDvh;   QMLhr= ratio*QMLhr
          endif

       enddo
       !---------------  End of iterative adjustment section.  ------------------!

       !========================================================================!
       !           Add all source/sink terms to all predicted moments:          !
       !========================================================================!

       ! Q-Source/Sink Terms:
       Q(i,k) = Q(i,k)  -QNUvi -QVDvi -QVDvs -QVDvg -QVDvh
       QC(i,k)= QC(i,k) -QCLci -QCLcs -QCLcg -QCLch -QFZci
       QR(i,k)= QR(i,k) -QCLri +QMLsr -QCLrs -QCLrg +QMLgr -QCLrh +QMLhr -QFZrh +QMLir
       QI(i,k)= QI(i,k) +QNUvi +QVDvi +QCLci -QCNig +QFZci -QCNis -QCLir -QCLis -QCLig   &
                        -QMLir -QCLih +QIMsi +QIMgi
       QG(i,k)= QG(i,k) +QCNsg +QVDvg +QCLcg -QCLgr-QMLgr -QCNgh -QIMgi +QCNig +QCLig    &
                        +Dirg*(QCLri+QCLir) +Dgrg*(QCLrg+QCLgr) +Dsrg*(QCLrs+QCLsr)
       QN(i,k)= QN(i,k) +QCNis +QVDvs +QCLcs -QCNsg -QMLsr -QIMsi -QCLsr +QCLis -QCLsh   &
                        +Dsrs*(QCLrs+QCLsr)
       QH(i,k)= QH(i,k) +Dirh*(QCLri+QCLir) -QMLhr +QVDvh +QCLch +Dsrh*(QCLrs+QCLsr)     &
                        +QCLih +QCLsh +QFZrh +QCLrh +QCNgh +Dgrh*(QCLrg+QCLgr)
       T(i,k)= T(i,k)   +LFP*(QCLci+QCLri+QCLcs+QCLrs+QFZci-QMLsr+QCLcg+QCLrg-QMLir      &
                        -QMLgr-QMLhr+QCLch+QCLrh+QFZrh) +LSP*(QNUvi+QVDvi+QVDvs+QVDvg    &
                        +QVDvh)

      if(QC(i,k)<epsQ)    then
         Q(i,k) = Q(i,k) + QC(i,k)
         T(i,k) = T(i,k) - LCP*QC(i,k)
         QC(i,k)= 0.
      endif
      if (QR(i,k)<epsQ)   then
         Q(i,k) = Q(i,k) + QR(i,k)
         T(i,k) = T(i,k) - LCP*QR(i,k)
         QR(i,k)= 0.
      endif
      if (QI(i,k)<epsQ)   then
         Q(i,k) = Q(i,k) + QI(i,k)
         T(i,k) = T(i,k) - LSP*QI(i,k)
         QI(i,k)= 0.
      endif
      if (QN(i,k)<epsQ)   then
         Q(i,k) = Q(i,k) + QN(i,k)
         T(i,k) = T(i,k) - LSP*QN(i,k)
         QN(i,k)= 0.
      endif
      if (QG(i,k)<epsQ)   then
         Q(i,k) = Q(i,k) + QG(i,k)
         T(i,k) = T(i,k) - LSP*QG(i,k)
         QG(i,k)= 0.
      endif
      if (QH(i,k)<epsQ)   then
         Q(i,k) = Q(i,k) + QH(i,k)
         T(i,k) = T(i,k) - LSP*QH(i,k)
         QH(i,k)= 0.
      endif
      Q(i,k)= max(Q(i,k), 0.)

      ENDIF  !if (activePoint)
    ENDDO
  ENDDO


  !----------------------------------------------------------------------------------!
  !                    End of ice phase microphysics (Part 2)                        !
  !----------------------------------------------------------------------------------!

  !----------------------------------------------------------------------------------!
  !                       PART 3: Warm Microphysics Processes                        !
  !                                                                                  !
  !  Equations for warm-rain coalescence based on Cohard and Pinty 2000a,b (QJRMS)   !
  !  Condensation/evaportaion equation based on Kong and Yau 1997 (Atmos-Ocean)      !
  !  Equations for rain reflectivity (ZR) based on Milbrandt and Yau 2005b (JAS)     !
  !----------------------------------------------------------------------------------!

  ! Warm-rain Coallescence:

 IF (warmphase_ON) THEN

  DO k= 2,nk
     DO i= 1,ni

        iDE   = 1./DE(i,k)
        RCAUTR= 0.;  CCACCR= 0.;  Dc= 0.;  iLAMc= 0.;  L  = 0.
        RCACCR= 0.;  CCSCOC= 0.;  Dr= 0.;  iLAMr= 0.;  TAU= 0.
        CCAUTR= 0.;  CRSCOR= 0.;  SIGc= 0.;  DrINIT= 0.
        iLAMc3= 0.;  iLAMc6= 0.;  iLAMr3= 0.;  iLAMr6= 0.

        rainPresent= (QRM(i,k)>epsQ)
        if (QCM(i,k)>epsQ) then
           iLAMc = ((DE(i,k)*(QCM(i,k)/N_c))*icexc9)**thrd
           iLAMc3= iLAMc*iLAMc*iLAMc
           iLAMc6= iLAMc3*iLAMc3
           Dc    = iLAMc*(GC2*iGC1)**thrd
           SIGc  = iLAMc*( GC3*iGC1- (GC2*iGC1)*(GC2*iGC1) )**sixth
           L     = 0.027*DE(i,k)*QCM(i,k)*(6.25e18*SIGc*SIGc*SIGc*Dc-0.4)
           if (SIGc>SIGcTHRS) TAU= 3.7/(DE(i,k)*(QCM(i,k))*(0.5e6*SIGc-7.5))
        endif

        if (rainPresent) then
           N_r= GR50*sqrt(sqrt(GR31/GR34*DE(i,k)*QRM(i,k)*icmr))
           Dr    = (DE(i,k)*(QRM(i,k)/N_r)*icmr)**thrd
           iLAMr = ((DE(i,k)*(QRM(i,k)/N_r))*icexr9)**thrd
           iLAMr3= iLAMr*iLAMr*iLAMr
           iLAMr6= iLAMr3*iLAMr3
        endif

        !  Autoconversion:
        if (QCM(i,k)>epsQ .and. SIGc>SIGcTHRS .and. autoconv_ON) then
           RCAUTR= min( max(L/TAU,0.), QCM(i,k)*idt )
           DrINIT= max(83.e-6, 12.6e-4/(0.5e6*SIGc-3.5))  !initiation regime Dr
           DrAUT = max(DrINIT, Dr)                     !init. or feeding DrAUT

           ! ---------------------------------------------------------------------------- !
           ! NOTE: The formulation for CCAUTR here (dNr/dt|initiation) does NOT follow
           !       eqn (18) in CP2000a, but rather it comes from the F90 code provided
           !       by J-P Pinty (subroutine: 'rain_c2r2.f90').
           !
           !       Similarly, the condition for the activation of accretion and self-
           !       collection depends on whether or not autoconversion is in the feeding
           !       regime.  This is apparent in the F90 code, but NOT in CP2000a.
           ! ---------------------------------------------------------------------------- !

        endif

        ! Accretion, rain self-collection, and collisional breakup:
        if (((QRM(i,k))>1.2*max(L,0.)*iDE.or.Dr>max(5.e-6,DrINIT)) .and. rainAccr_ON     &
             .and. rainPresent) then

           !  Accretion:                                                      !{CP00a eqn(22)}
           if (QCM(i,k)>epsQ.and.L>0.) then
              if (Dr.ge.100.e-6) then
                 RCACCR = cmr*iDE*KK1*(N_c*N_r)*iLAMc3*(GC3*iGC1*iLAMc3+GC2*iGC1*GR34*   &
                          iGR31*iLAMr3)
              else
!++  The following calculation of RCACCR avoids overflow:
                 tmp1   = cmr*iDE
                 tmp2   = KK2*(N_c*N_r)*iLAMc3
                 RCACCR = tmp1 * tmp2
                 tmp1   = GC4*iGR31
                 tmp1   = (tmp1)*iLAMc6
                 tmp2   = GC2*iGC1
                 tmp2   = tmp2*GR37*iGR31
                 tmp2   = (tmp2)*iLAMr6
                 RCACCR = RCACCR * (tmp1 + tmp2)
!++
              endif
              RCACCR = min(RCACCR,(QC(i,k))*idt)
            endif

        endif  !accretion/self-collection/breakup

        ! Prevent overdepletion of cloud:
        source= (QC(i,k))
        sink  = (RCAUTR+RCACCR)*dt
        if (sink>source) then
           ratio = source/sink
           RCAUTR= ratio*RCAUTR
           RCACCR= ratio*RCACCR
        endif

        ! Apply tendencies:
        QC(i,k)= max(0., QC(i,k)+(-RCAUTR-RCACCR)*dt )
        QR(i,k)= max(0., QR(i,k)+( RCAUTR+RCACCR)*dt )

     ENDDO
  ENDDO

  ! Condensation/Evaporation:

  DO k=1,nk
     DO i=1,ni

        iDE     = 1./DE(i,k)
        DEo     = DE(i,nk)
        gam     = sqrt(DEo*iDE)
        QSS(i,k)= sngl(FOQSA(T(i,k), PS(i)*S(i,k)))  ! Re-calculates QS with new T (w.r.t. liquid)
       !----
       !The following removes a fraction of the supersaturation water vapor.  The purpose of this
       ! is to reduce the amount of condensed water, ultimately to reduct the total precipitation.
       ! Note -- there is no physical justification for this within the context of the microphysics.
       ! This adjustment is done entirely to reduce the precipitation bias in the GEM-LAM-2.5.
        Q(i,k)= Q(i,k) - max(0., qReducFact*(Q(i,k)-QSS(i,k)))
       !----
        ssat    = Q(i,k)/QSS(i,k)-1.
        Tc      = T(i,k)-TRPL
        Cdiff   = max(1.62e-5, (2.2157e-5 + 0.0155e-5*Tc)) *1.e5/(S(i,k)*PS(i))
        MUdyn   = max(1.51e-5, (1.7153e-5 + 0.0050e-5*Tc))
        MUkin   = MUdyn*iDE
        iMUkin  = 1./MUkin
        Ka      = max(2.07e-2, (2.3971e-2 + 0.0078e-2*Tc))
        ScTHRD  = (MUkin/Cdiff)**thrd ! i.e. Sc^(1/3)

        !Condensation/evaporation:
        ! Capacity of evap/cond in one time step is determined by saturation
        ! adjustment technique [KY97 App.A].  Equation for rain evaporation rate
        ! comes from CP00a.  Explicit condensation rate is not considered
        ! (as it is in Z85), but rather complete removal of supersaturation
        ! is assumed.

        X= Q(i,k)-QSS(i,k)
        rainPresent= (QR(i,k)>epsQ)
        IF(X>0. .or. QC(i,k)>epsQ .or. rainPresent) THEN
           tmp1 = T(i,k)-35.86
           X    = X/(1.+ck5*QSS(i,k)/(tmp1*tmp1))
           if (X<(-QC(i,k))) then
              D= 0.
              if(rainPresent) then

                 if(QM(i,k)<QSW(i,k)) then
                    MUkin = (1.715e-5+5.e-8*Tc)*iDE
                    iMUkin= 1./MUkin
                    ! Rain evap: (F94)
                    iLAMr  = sqrt(sqrt( (QR(i,k)*DE(i,k))/(GR34*cmr*No_r) ))
                    VENTr= Avx*GR32*iLAMr*iLAMr + Bvx*ScTHRD*sqrt(gam*afr*iMUkin)*GR17*  &
                           iLAMr**cexr6                         !optimized for alpha_r=0
                      !note: There is an error in the documented MY05a_eq(8) for VENTx
                      !      (The above code is correct)
                    ABw  = CHLF*CHLF/(Ka*RGASV*T(i,k)*T(i,k))+1./(DE(i,k)*(QSS(i,k))*    &
                           Cdiff)
                    QREVP= -dt*(PI2*ssat*No_r*VENTr/ABw)
                    if ((QR(i,k))>QREVP) then             !Note: QREVP is [(dQ/dt)*dt]
                       DEL= -QREVP
                    else
                       DEL= -QR(i,k)
                    endif
                    D= max(X+QC(i,k), DEL)
                 endif  !QM< QSM
              endif   !QR<eps
              X= D - QC(i,k)

              QR(i,k)= QR(i,k) + D
              ! The above expression of (dNr/dt)|evap is from Ferrier, 1994.
              ! In CP2000a, Nr is not affected by evap. (except if Qr goes to zero).
              QC(i,k)= 0.
              T(i,k) = T(i,k) + LCP*X
              Q(i,k) = Q(i,k) - X

           else  ![if(X >= -QC)]

              ! All supersaturation is removed (condensed onto cloud field).
              T(i,k)  = T(i,k)  + LCP*X
              Q(i,k)  = Q(i,k)  - X
              QC(i,k) = QC(i,k) + X

              ! Homogeneous freezing of cloud to ice:
              !  Note:  This needs to be calculated here, as well as in the ice-phase
              !         section, in case water condenses at a very low temperature.
            ! Simplified:
              if (QC(i,k)>epsQ .and. T(i,k)<238.15 .and. icephase_ON) then
                  FRAC    = 1.  !if T<-35C
                  QFZci   = FRAC*QC(i,k)
                  QI(i,k) = QI(i,k) + QFZci
                  QC(i,k) = QC(i,k) - QFZci
                  T(i,k)  = T(i,k)  + QFZci*LFP
              endif !homogeneous freezing

           endif

        ENDIF

       !Protect against negative values due to overdepletion:
        if (QR(i,k)<epsQ)  then
           Q(i,k) = Q(i,k) + QR(i,k)
           T(i,k) = T(i,k) - QR(i,k)*LCP
           QR(i,k)= 0.
        endif

     ENDDO
  ENDDO    !cond/evap [k-loop]

 ENDIF  !if warmphase_ON

  !----------------------------------------------------------------------------------!
  !                    End of warm-phase microphysics (Part 3)                       !
  !----------------------------------------------------------------------------------!

  !----------------------------------------------------------------------------------!
  !                            PART 4:  Sedimentation                                !
  !----------------------------------------------------------------------------------!

  !----------------------------------------------------------------------------------!
  ! Sedimentation is computed using a modified version of the box-Lagrangian         !
  ! scheme (blg3sedi.ftn).  The modifications (from blg2.ftn) are hard-wired         !
  ! for use with multimom_xx.ftn90.  Sedimentation is only computed for columns      !
  ! containing non-zero hydrometeor quantities (at at least one level).              !
  !----------------------------------------------------------------------------------!


 IF (sedi_ON) THEN
!  call tmg_start0(51,'mmSEDI')

   RT_rn1= 0.;  RT_rn2= 0.;  RT_fr1= 0.;  RT_fr2= 0.;  RT_sn1= 0.
   RT_sn2= 0.;  RT_sn3= 0.;  RT_pe1= 0.;  RT_pe2= 0.;  RT_peL= 0.

!--  RAIN sedimentation:
!  call tmg_start0(91,'mmSEDI_r')
   call SEDI_main_1(QR,1,T,DE,gamfact_r,epsQr_sedi,afr,bfr,icmr,dmr,dtr,cr6,ckQr1,ckQr2, &
                    icexr9,npassr,ni,nk,VrMax,DrMax,DZ,RT_rn1,No_r,ktop_sedi,            &
                    massFlux3D=massFlux3D_r)
!  call tmg_stop0(91)


!--  ICE  sedimentation:
!  call tmg_start0(92,'mmSEDI_i')
   call SEDI_main_1(QI,2,T,DE,gamfact,epsQi_sedi,afi,bfi,icmi,dmi,dti,ci6,ckQi1,ckQi2,   &
                   ckQi4,npassi,ni,nk,ViMax,DiMax,DZ,RT_sn1,-99.,ktop_sedi,              &
                   massFlux3D=massFlux3D_i)
!  call tmg_stop0(92)


!--  SNOW sedimentation:
!  call tmg_start0(93,'mmSEDI_s')
   call SEDI_main_1(QN,3,T,DE,gamfact,epsQs_sedi,afs,bfs,icms,dms,dts,cs6,ckQs1,ckQs2,   &
                    ckQs4,npasss,ni,nk,VsMax,DsMax,DZ,RT_sn2,-99.,ktop_sedi,             &
                    massFlux3D=massFlux3D_s)
!  call tmg_stop0(93)


!--  GRAUPEL sedimentation:
!  call tmg_start0(94,'mmSEDI_g')
   call SEDI_main_1(QG,4,T,DE,gamfact,epsQg_sedi,afg,bfg,icmg,dmg,dtg,cg6,ckQg1,ckQg2,   &
                    ckQg4,npassg,ni,nk,VgMax,DgMax,DZ,RT_sn3,No_g,ktop_sedi,             &
                    massFlux3D=massFlux3D_g)
!  call tmg_stop0(94)


!--  HAIL sedimentation:
!  call tmg_start0(95,'mmSEDI_h')
   call SEDI_main_1(QH,5,T,DE,gamfact,epsQh_sedi,afh,bfh,icmh,dmh,dth,ch6,ckQh1,ckQh2,   &
                    ckQh4,npassh,ni,nk,VhMax,DhMax,DZ,RT_pe1,No_h,ktop_sedi)
!  call tmg_stop0(95)

!---  End of sedimentation for each category --------!

 ! Compute MASS fluxes at surface:  (liquid-equivalent precipitation rates) [kg/m^2/s]
   RT_rn1 = RT_rn1 *ck7
   RT_sn1 = RT_sn1 *ck7
   RT_sn2 = RT_sn2 *ck7
   RT_sn3 = RT_sn3 *ck7
   RT_pe1 = RT_pe1 *ck7

 ! Compute VOLUME fluxes at surface:  (liquid-equivalent precipitation rates) [m/s]
   RT_rn1 = RT_rn1 *idew
   RT_sn1 = RT_sn1 *idew
   RT_sn2 = RT_sn2 *idew
   RT_sn3 = RT_sn3 *idew
   RT_pe1 = RT_pe1 *idew

  !Compute sum of solid (unmelted) volume fluxes [m3 (bulk hydrometeor) m-2 (sfc area) s-1]:
  !  (this is the precipitation rate for UNmelted total "snow" [i+s+g])
  !        In 'calcdiagn.ftn', the total solid (excluding hail) precipitation, SN, is computed
  !        from the sum of the liq-eq.vol fluxes, tss_sn1 + tss_sn2 + tss_sn3.  With the
  !        accumulation of sn_bulk (in 'calcdiag.ftn'), the solid-liquid ratio for the total
  !        snow (i+s+g) can be compute as SNblk/SN.

   RT_snd = massFlux3D_i(:,nk)/dei + massFlux3D_g(:,nk)/deg
   do i= 1,ni
      if (QN(i,nk)>epsQ) then
         if (T(i,nk)<TRPL) then
            des_mlt = des
         else
           !estimate bulk density of melting snow by approximating liquid fraction by QR/(QR+QN):
            des_mlt = (QN(i,nk)*des + QR(i,nk)*dew)/(QN(i,nk)+QR(i,nk))
         endif
         RT_snd(i) = RT_snd(i) + massFlux3D_s(i,nk)/des_mlt
      endif
   enddo


 !++++
 ! Diagnose surface precipitation types:

 ! The following involves diagnostic conditions to determine surface precipitation rates
 ! for various precipitation elements defined in Canadian Meteorological Operational Internship
 ! Program module 3.1 (plus one for large hail) based on the sedimentation rates of the five
 ! sedimenting hydrometeor categories.
 !
 ! With the diagnostics shut off (precipDiag_ON=.false.), 5 rates are just the 5 category
 ! rates, with the other 6 rates just 0.  The model output variables will have:
 !
 !   total liquid = RT_rn1 [RAIN]
 !   total solid  = RT_sn1 [ICE] + RT_sn2 [SNOW] + RT_sn3 [GRAUPEL] + RT_pe1 [HAIL]
 !
 ! With the diagnostics on, the 5 sedimentation rates are partitioned into 9 rates,
 ! with the following model output variable:
 !
 !  total liquid = RT_rn1 [liquid rain] + RT_rn2 [liquid drizzle]
 !  total solid  = RT_fr1 [freezing rain] + RT_fr2 [freezing drizzle] + RT_sn1 [ice crystals] +
 !                 RT_sn2 [snow] + RT_sn3 [graupel] + RT_pe1 [ice pellets] + RT_pe2 [hail]
 !
 ! NOTE:  - The above total liquid/solid rates are computed in 'calcdiag.ftn' (as R2/R4).
 !        - RT_peL [large hail] is a sub-set of RT_pe2 [hail] and is thus not added to the total
 !          solid precipitation.

   IF (precipDiag_ON) THEN
      DO i= 1,ni

         DE(i,nk)= S(i,nk)*PS(i)/(RGASD*T(i,nk))

      !rain vs. drizzle:
         N_r= (No_r*GR31)**(3./(4.+alpha_r))*(GR31/GR34*DE(i,nk)*QR(i,nk)*icmr)**        &
               ((1.+alpha_r)/(4.+alpha_r))             !i.e. NR = f(No_r,QR)
         if (QR(i,nk)>epsQ .and. N_r>epsN) Dm_r(i,nk)= (DE(i,nk)*icmr*QR(i,nk)/N_r)**thrd
         if (Dm_r(i,nk)>Dr_large) then  !size threshold to distinguish rain/drizzle
            RT_rn2(i)= RT_rn1(i);   RT_rn1(i)= 0.
         endif

      !liquid vs. freezing rain or drizzle:
         if (T(i,nk)<TRPL) then
            RT_fr1(i)= RT_rn1(i);   RT_rn1(i)= 0.
            RT_fr2(i)= RT_rn2(i);   RT_rn2(i)= 0.
         endif

      !ice pellets vs. hail:
         if (T(i,nk)>(TRPL+5.0)) then
         !note: The condition (T_sfc<5C) for ice pellets is a simply proxy for the presence
         !      of a warm layer aloft, though which falling snow or graupel will melt to rain,
         !      over a sub-freezing layer, where the rain will freeze into the 'hail' category
            RT_pe2(i)= RT_pe1(i);   RT_pe1(i)= 0.
         endif

      !large hail:
         if (QH(i,nk)>epsQ) then
            N_h= (No_h*GH31)**(3./(4.+alpha_h))*(GH31*iGH34*DE(i,nk)*QH(i,nk)*           &
                 icmh)**((1.+alpha_h)/(4.+alpha_h))   !i.e. Nh = f(No_h,Qh)
            Dm_h(i,nk)= (DE(i,nk)*icmh*QH(i,nk)/N_h)**thrd
            if (DM_h(i,nk)>Dh_large) RT_peL(i)= RT_pe2(i)
            !note: large hail (RT_peL) is a subset of the total hail (RT_pe2)
         endif

      ENDDO
   ENDIF    ! if (precipDiag_ON)
!++++

!  call tmg_stop0(51)
 ENDIF  ! if (sedi_ON)

 where (Q<0.) Q= 0.

 !-----------------------------------------------------------------------------------!
 !                     End of sedimentation calculations (Part 4)                    !
 !-----------------------------------------------------------------------------------!


 !===================================================================================!
 !                             End of microphysics scheme                            !
 !===================================================================================!

 !-----------------------------------------------------------------------------------!
 !                    Calculation of diagnostic output variables:                    !

  IF (calcDiag) THEN

   !For reflectivity calculations:
     cxr =            icmr*icmr   !for rain
     cxi = 1./fdielec*icmr*icmr   !for all frozen categories
     Gzr = (6.+alpha_r)*(5.+alpha_r)*(4.+alpha_r)/((3.+alpha_r)*(2.+alpha_r)*(1.+alpha_r))
     Gzi = (6.+alpha_i)*(5.+alpha_i)*(4.+alpha_i)/((3.+alpha_i)*(2.+alpha_i)*(1.+alpha_i))
     Gzs = (6.+alpha_s)*(5.+alpha_s)*(4.+alpha_s)/((3.+alpha_s)*(2.+alpha_s)*(1.+alpha_s))
     Gzg = (6.+alpha_g)*(5.+alpha_g)*(4.+alpha_g)/((3.+alpha_g)*(2.+alpha_g)*(1.+alpha_g))
     Gzh = (6.+alpha_h)*(5.+alpha_h)*(4.+alpha_h)/((3.+alpha_h)*(2.+alpha_h)*(1.+alpha_h))

     do k= 1,nk
       do i= 1,ni
          DE(i,k)= S(i,k)*PS(i)/(RGASD*T(i,k))
          tmp9= DE(i,k)*DE(i,k)

        !Compute Nt_x (single-moment)
          N_r= (No_r*GR31)**(3./(4.+alpha_r))*(GR31/GR34*DE(i,k)*     &
                    QR(i,k)*icmr)**((1.+alpha_r)/(4.+alpha_r))      !i.e. NRM = f(No,QRM)
          N_i= 5.*exp(0.304*(TRPL-max(233.,T(i,k))))                !Cooper eqn.
          No_s= min(2.e+8, 2.e+6*exp(-0.12*min(-0.001,T(i,k)-TRPL)))!Thompson et al. (2004)
          N_s= (No_s*GS31)**(3./(4.+alpha_s))*(GS31/GS34*DE(i,k)*     &
                    QN(i,k)*icms)**((1.+alpha_s)/(4.+alpha_s))      !i.e. NXM = f(No,QXM)
          N_g= (No_g*GG31)**(3./(4.+alpha_g))*(GG31/GG34*DE(i,k)*     &
                    QG(i,k)*icmg)**((1.+alpha_g)/(4.+alpha_g))      !i.e. NXM = f(No,QXM)
          N_h= (No_h*GH31)**(3./(4.+alpha_h))*(GH31*iGH34*DE(i,k)*    &
                    QH(i,k)*icmh)**((1.+alpha_h)/(4.+alpha_h))      !i.e. NXM = f(No,QXM)

        !Total equivalent reflectivity:     (units of [dBZ])
           tmp1= 0.;  tmp2= 0.;  tmp3= 0.;  tmp4= 0.;  tmp5= 0.
           if (QR(i,k)>epsQ .and. N_r>epsN) tmp1 = cxr*Gzr*tmp9*QR(i,k)*QR(i,k)/N_r
           if (QI(i,k)>epsQ .and. N_i>epsN) tmp2 = cxi*Gzi*tmp9*QI(i,k)*QI(i,k)/N_i
           if (QN(i,k)>epsQ .and. N_s>epsN) tmp3 = cxi*Gzs*tmp9*QN(i,k)*QN(i,k)/N_s
           if (QG(i,k)>epsQ .and. N_g>epsN) tmp4 = cxi*Gzg*tmp9*QG(i,k)*QG(i,k)/N_g
           if (QH(i,k)>epsQ .and. N_h>epsN) tmp5 = cxi*Gzh*tmp9*QH(i,k)*QH(i,k)/N_h
          !Modifiy dielectric constant for melting ice-phase categories:
           if ( T(i,k)>TRPL) then
             tmp2= tmp2*fdielec
             tmp3= tmp3*fdielec
             tmp4= tmp4*fdielec
             tmp5= tmp5*fdielec
           endif
           ZET(i,k) = tmp1 + tmp2 + tmp3 + tmp4 + tmp5   != Zr+Zi+Zs+Zg+Zh
           if (ZET(i,k)>0.) then
              ZET(i,k)= 10.*log10((ZET(i,k)*Zfact))      !convert to dBZ
           else
              ZET(i,k)= minZET
           endif
           ZET(i,k)= max(ZET(i,k),minZET)
           ZEC(i)= max(ZEC(i),ZET(i,k))  !composite (max in column) of ZET

         !Mean-mass diameters:  (units of [m])
           Dm_c(i,k)= 0.;   Dm_r(i,k)= 0.;   Dm_i(i,k)= 0.
           Dm_s(i,k)= 0.;   Dm_g(i,k)= 0.;   Dm_h(i,k)= 0.
           if(QC(i,k)>epsQ .and. N_c>epsN) Dm_c(i,k)= (DE(i,k)*icmr*QC(i,k)/N_c)**thrd
           if(QR(i,k)>epsQ .and. N_r>epsN) Dm_r(i,k)= (DE(i,k)*icmr*QR(i,k)/N_r)**thrd
           if(QI(i,k)>epsQ .and. N_i>epsN) Dm_i(i,k)= (DE(i,k)*icmi*QI(i,k)/N_i)**thrd
           if(QN(i,k)>epsQ .and. N_s>epsN) Dm_s(i,k)= (DE(i,k)*icms*QN(i,k)/N_s)**thrd
           if(QG(i,k)>epsQ .and. N_g>epsN) Dm_g(i,k)= (DE(i,k)*icmg*QG(i,k)/N_g)**thrd
           if(QH(i,k)>epsQ .and. N_h>epsN) Dm_h(i,k)= (DE(i,k)*icmh*QH(i,k)/N_h)**thrd

         !Supercooled liquid water:
           SLW(i,k)= 0.
           if (T(i,k)<TRPL) SLW(i,k)= DE(i,k)*(QC(i,k)+QR(i,k))   !(units of [kg m-3])

         !Visibility:

          !VIS1:  component through liquid cloud (fog) [m]
          ! (following parameterization of Gultepe and Milbrandt, 2007)
           tmp1= QC(i,k)*DE(i,k)*1.e+3             !LWC [g m-3]
           tmp2= N_c*1.e-6                         !Nc  [cm-3]
           if (tmp1>0.005 .and. tmp2>1.) then
              VIS1(i,k)= 1000.*(1.13*(tmp1*tmp2)**(-0.51)) !based on FRAM [GM2007, eqn (4)
             !VIS1(i,k)= min(maxVIS, (tmp1*tmp2)**(-0.65)) !based on RACE [GM2007, eqn (3)
           else
              VIS1(i,k)= 3.*maxVIS  !gets set to maxVIS after calc. of VIS
           endif

          !VIS2: component through rain  !based on Gultepe and Milbrandt, 2008, Table 2 eqn (1)
           tmp1= massFlux3D_r(i,k)*idew*3.6e+6                        !rain rate [mm h-1]
           if (tmp1>0.01) then
              VIS2(i,k)= 1000.*(-4.12*tmp1**0.176+9.01)               ![m]
           else
              VIS2(i,k)= 3.*maxVIS
           endif

          !VIS3: component through snow  !based on Gultepe and Milbrandt, 2008, Table 2 eqn (6)
           tmp1= massFlux3D_s(i,k)*idew*3.6e+6                        !snow rate, liq-eq [mm h-1]
           if (tmp1>0.01) then
              VIS3(i,k)= 1000.*(1.10*tmp1**(-0.701))                  ![m]
           else
              VIS3(i,k)= 3.*maxVIS
           endif

          !VIS:  visibility due to reduction from all components 1, 2, and 3
          !      (based on sum of extinction coefficients and Koschmieders's Law)
           VIS(i,k) = min(maxVIS, 1./(1./VIS1(i,k) + 1./VIS2(i,k) + 1./VIS3(i,k)))
           VIS1(i,k)= min(maxVIS, VIS1(i,k))
           VIS2(i,k)= min(maxVIS, VIS2(i,k))
           VIS3(i,k)= min(maxVIS, VIS3(i,k))

!VIS(i,k)= massFlux3D_r(i,k)

        enddo  !i-loop
     enddo     !k-loop

    !Diagnostic levels:
     h_CB = noVal_h_XX   !height (AGL) of cloud base
     h_SN = noVal_h_XX   !height (AGL) of snow level [conventional snow (not just QN>0.)]
     h_ML1= noVal_h_XX   !height (AGL) of melting level [first 0C isotherm from ground]
     h_ML2= noVal_h_XX   !height (AGL) of melting level [first 0C isotherm from top]
                           !note: h_ML2 = h_ML1 implies only 1 melting level
     tmp1= 1./GRAV
     do i= 1,ni
        CB_found= .false.;   SN_found= .false.;   ML_found= .false.
        do k= nk,2,-1
          !cloud base:
           if ((QC(i,k)>epsQ.or.QI(i,k)>epsQ) .and. .not.CB_found) then
              h_CB(i) = GZ(i,k)*tmp1
              CB_found= .true.
           endif
          !snow level:
           if ( ((QN(i,k)>epsQ .and. Dm_s(i,k)>minSnowSize)  .or.                        &
                 (QG(i,k)>epsQ .and. Dm_g(i,k)>minSnowSize)) .and. .not.SN_found) then
              h_SN(i) = GZ(i,k)*tmp1
              SN_found= .true.
           endif
          !melting level: (height of lowest 0C isotherm)
           if (T(i,k)>TRPL .and. T(i,k-1)<TRPL .and. .not.ML_found) then
              h_ML1(i) = GZ(i,k)*tmp1
              ML_found= .true.
           endif
        enddo
     enddo
     do i= 1,ni
        ML_found= .false.  !from top
        do k= 2,nk
          !melting level: (height of highest 0C isotherm)
           if (T(i,k)>TRPL .and. T(i,k-1)<TRPL .and. .not.ML_found) then
              h_ML2(i) = GZ(i,k)*tmp1
              ML_found= .true.
           endif
        enddo
     enddo

  ENDIF
 !                                                                                   !
 !-----------------------------------------------------------------------------------!

 !-----------------------------------------------------------------------------------!
 !   Compute the tendencies of  T, Q, QC, etc. (to be passed back to model dynamics) !
 !   and reset the fields to their initial (saved) values at time {*}:               !

      do k= 1,nk
         do i= 1,ni
            tmp1=T_TEND(i,k);  T_TEND(i,k)=(T(i,k) -T_TEND(i,k))*idt;  T(i,k) = tmp1
            tmp1=Q_TEND(i,k);  Q_TEND(i,k)=(Q(i,k) -Q_TEND(i,k))*idt;  Q(i,k) = tmp1
            tmp1=QCTEND(i,k);  QCTEND(i,k)=(QC(i,k)-QCTEND(i,k))*idt;  QC(i,k)= tmp1
            tmp1=QRTEND(i,k);  QRTEND(i,k)=(QR(i,k)-QRTEND(i,k))*idt;  QR(i,k)= tmp1
            tmp1=QITEND(i,k);  QITEND(i,k)=(QI(i,k)-QITEND(i,k))*idt;  QI(i,k)= tmp1
            tmp1=QNTEND(i,k);  QNTEND(i,k)=(QN(i,k)-QNTEND(i,k))*idt;  QN(i,k)= tmp1
            tmp1=QGTEND(i,k);  QGTEND(i,k)=(QG(i,k)-QGTEND(i,k))*idt;  QG(i,k)= tmp1
            tmp1=QHTEND(i,k);  QHTEND(i,k)=(QH(i,k)-QHTEND(i,k))*idt;  QH(i,k)= tmp1
         enddo
      enddo
 !                                                                                   !
 !-----------------------------------------------------------------------------------!

 END SUBROUTINE MYSMOM_MAIN
!===================================================================================================!

END MODULE my_smom_mod