updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / chem / module_cam_mam_wetscav.F
blob4ac5fc056d4f2eeb6501fc8feb894c157afbec83
1 !Updated to CESM1.0.3 (CAM5.1.01) by Balwinder.Singh@pnnl.gov
2 #define MODAL_AERO
3 module module_cam_mam_wetscav
4   !==================================================================================
5   !Future Modifications:
6   !---------------------
7   !1. CME (or QME:Net condensation rate) is NOT used in this subroutine, although it 
8   !   is being passed to the wet scavenging parameterizations. The value stored in
9   !   QME3D variable is not correct as it lacks the contribution of CMELIQ from CAM's 
10   !   Macrophysics
11   !2. CAM_OUT datastructure of CAM provides information to the surface model of CAM. 
12   !   It has NOT been implemented in WRF.
13   !3. LCHNK is being passed to the parameterization and it is used by CAM's OUTFLD 
14   !   calls. OUTFLD calls are just dummy calls for WRF. Therefore LCHNK is not used 
15   !   for any meaningful calculation which can alter the answers from the 
16   !   parameterization
17   !4. Currently the CAM's wet scavenging is forced to ONLY work with the CAMMGMP 
18   !   microphysics
19   !5. LAT and LON variables are ONLY used in chemistry calls in CAM, which is NOT 
20   !   implemented in WRF yet
21   !6. Variable "num_mz_aerosols" in module_cam_mam_mz_aerosols_intr.F parameterization 
22   !   Fortran file is HARDWIRED to be equal to 'pcnst' (basically, it has to be greater 
23   !   than zero for the parameterization to run)
24   !
25   !==================================================================================
26   !
27   !!---------------------------------------------------------------------------------
28   !! Module to interface the aerosol parameterizations with CAM
29   !! Phil Rasch, Jan 2003
30   !  Ported to WRF by Balwinder.Singh@pnnl.gov
31   !!---------------------------------------------------------------------------------
32   !
33   use shr_kind_mod, only: r8 => shr_kind_r8, cl => shr_kind_cl
34   
35   implicit none
36   private
37   save
38   !
39   !! Public interfaces
40   !
41   
42   public :: wetscav_cam_mam_driver_init ! initialization subroutine
43   public :: wetscav_cam_mam_driver      ! interface to wet deposition 
44   
45   integer      :: itf, jtf, ktf, pverp
46   real(r8)     :: dp1, dp2
47   
48   !===============================================================================
49 contains
50   !===============================================================================
51   
52   subroutine wetscav_cam_mam_driver(itimestep,p_hyd,p8w,t_phy,             &
53        dgnum4d,dgnumwet4d,dlf3d,dlf2_3d,dtstep,qme3d,                      &
54        prain3d,nevapr3d,rate1ord_cw2pr_st3d,shfrc3d,cmfmc3d,cmfmc2_3d,     &
55        evapcsh3d,icwmrsh3d,rprdsh3d,evapcdp3d,icwmrdp3d,rprddp3d,          &
56        qs_curr, f_ice_phy,f_rain_phy,config_flags,cldfra_mp_all,cldfrai,   &
57        cldfral,cldfra,is_CAMMGMP_used,                                     &
58        ids,ide, jds,jde, kds,kde,                                          &
59        ims,ime, jms,jme, kms,kme,                                          &
60        its,ite, jts,jte, kts,kte,                                          &
61        !intent-inout
62        qv_curr,qc_curr,qi_curr,ni3d,nc3d,chem,                             &
63        !intent-out             
64        fracis3D                                                            )
65     
66     !!----------------------------------------------------------------------- 
67     !! 
68     !! Purpose: 
69     !! Interface to wet processing of all aerosols
70     !! 
71     !! Method: 
72     !!  use a modified version of the scavenging parameterization described in
73     !!     Barth et al, 2000, JGR (sulfur cycle paper)
74     !!     Rasch et al, 2001, JGR (INDOEX paper)
75     !! 
76     !! Author: Phil Rasch
77     !  Ported to WRF by Balwinder.Singh@pnnl.gov
78     !! 
79     !!-----------------------------------------------------------------------
80     use module_mp_cammgmp_driver,  only: physics_update, physics_ptend_init
81     use module_cam_support,        only: pcnst =>pcnst_runtime, pcols, pver
82     use module_state_description,  only: num_chem, param_first_scalar,F_QC, F_QI, F_QV, F_QS, &
83          CAMZMSCHEME, CAMUWSHCUSCHEME
84     use module_data_cam_mam_asect, only: lptr_chem_to_q, lptr_chem_to_qqcw, factconv_chem_to_q, waterptr_aer
85     use wetdep,                    only: clddiag
86 #if ( defined TROPCHEM || defined MODAL_AERO )
87     use mz_aerosols_intr,          only: mz_aero_wet_intr
88 #endif
89 #if ( defined MODAL_AERO )
90     use modal_aero_data,           only: ntot_amode
91 #endif
92     use module_configure,          only: grid_config_rec_type
93     use infnan,                    only: nan
94     !-----------------------------------------------------------------------
95     implicit none
96     !-----------------------------------------------------------------------
97     !
98     ! Arguments:
99     !
100     !
101     TYPE(grid_config_rec_type),  INTENT(IN   )    :: config_flags
103     logical, intent(in) :: is_CAMMGMP_used
105     integer, intent(in) :: itimestep                           !Time step number
106     integer, intent(in) :: ids,ide, jds,jde, kds,kde
107     integer, intent(in) :: ims,ime, jms,jme, kms,kme
108     integer, intent(in) :: its,ite, jts,jte, kts,kte
110     real, intent(in)    :: dtstep                              !Time step in seconds(s)
111     
112     !3d real arrays
113     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: dlf3d          !Detraining cloud water tendency from convection
114     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: dlf2_3d        !dq/dt due to export of cloud water into environment by shallow convection(kg/kg/s)
115     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: p8w            !Hydrostatic Pressure at level interface (Pa)
116     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: p_hyd          !Hydrostatic pressure(Pa)
117     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: t_phy          !Temperature (K)
118     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: qme3d          !Net condensation rate (kg/kg/s) *NOTE: Doesn't include contribution from CMELIQ of MACRO-PHYSICS*
119     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: prain3d        !Rate of conversion of condensate to precipitation (kg/kg/s)
120     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: nevapr3d       !Evaporation rate of rain + snow (kg/kg/s)
121     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: rate1ord_cw2pr_st3d !1st order rate for direct conversion of strat. cloud water to precip (1/s)
122     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: shfrc3d        !Shallow cloud fraction
123     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: cmfmc3d        !Deep + Shallow Convective mass flux [ kg /s/m^2 ]
124     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: cmfmc2_3d      !Shallow convective mass flux [ kg/s/m^2 ]
125     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: evapcsh3d      !Evaporation of shallow convection precipitation (kg/kg/s)
126     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: icwmrsh3d      !Shallow cumulus in-cloud water mixing ratio (kg/m2)
127     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: rprdsh3d       !dq/dt due to deep and shallow convective rainout(kg/kg/s)
128     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: evapcdp3d      !Evaporation of deep convective precipitation (kg/kg/s)
129     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: icwmrdp3d      !Deep Convection in-cloud water mixing ratio (kg/m2)
130     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: rprddp3d       !dq/dt due to deep convective rainout (kg/kg/s)
131     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: qs_curr        !Snow mixing ratio         -Cloud ice  (kg/kg)
132     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: F_ICE_PHY      !Fraction of ice
133     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: F_RAIN_PHY     !Fraction of rain
134     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: cldfra_mp_all
135     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: cldfrai
136     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: cldfral
137     real, dimension( ims:ime, kms:kme, jms:jme ), intent(in) :: cldfra
138     
139     !4D-Intent-in-real array
140     real, dimension( ims:ime, kms:kme, jms:jme, ntot_amode ), intent(in) :: dgnum4d    !4-dimensional Number mode diameters
141     real, dimension( ims:ime, kms:kme, jms:jme, ntot_amode ), intent(in) :: dgnumwet4d !4-dimensional Number mode diameters
143     !3D-Intent-inout-real array
144     real, dimension( ims:ime, kms:kme, jms:jme ), intent(inout) :: qv_curr     !Water vapor mixing ratio - Moisture  (kg/kg)
145     real, dimension( ims:ime, kms:kme, jms:jme ), intent(inout) :: qc_curr     !Cloud water mixing ratio - Cloud liq (kg/kg)
146     real, dimension( ims:ime, kms:kme, jms:jme ), intent(inout) :: qi_curr     !Ice mixing ratio         - Cloud ice  (kg/kg)
147     real, dimension( ims:ime, kms:kme, jms:jme ), intent(inout) :: ni3d        !Cloud ice number concentration (#/kg)
148     real, dimension( ims:ime, kms:kme, jms:jme ), intent(inout) :: nc3d        !Cloud water  number concentration (#/kg)
151     !4D-Intent-inout-real array
152     real, dimension( ims:ime, kms:kme, jms:jme, num_chem ),   intent(inout) :: chem     !Chem array
154     !3D-Intent-out-real array
155     real, dimension( ims:ime, kms:kme, jms:jme, pcnst ), intent(out) :: fracis3d      !Fraction of transported species that are insoluble
157     !Local variable specific to CAM 
158     real(r8) :: dt                      !Time step
159     integer  :: nstep                   !Time step number      
160     real(r8) :: clds(pcols,kte)         !Stratiform cloud fraction
161     
162 #if ( defined MODAL_AERO )
163     integer  i, k
164 #endif
165     
166     character*24 :: ptend_name            !ptend%name in CAM5 - used in parameterization
167     
168     logical      :: ptend_ls              !ptend%ls   in CAM5 - used for calling physics_update subroutine
169     logical      :: ptend_lq(pcnst)       !ptend%lq   in CAM5
170     
171     integer      :: iw, kw, jw, itsm1
172     integer      :: itile_len, ktep1
173     integer      :: kflip, l, l2, p1st
174     integer      :: imode, kcam
175     integer      :: ncol                  !state%ncol
176     integer      :: lchnk                 !state%lchnk
177     
178     real(r8)     :: dp, multFrc
179     
180     real(r8)     :: cldc(pcols,kte)           ! convective cloud fraction
181     real(r8)     :: cldv(pcols,kte)           ! cloudy volume undergoing scavenging
182     real(r8)     :: cldvcu(pcols,kte)         ! Convective precipitation area at the top interface of current layer
183     real(r8)     :: cldvst(pcols,kte)         ! Stratiform precipitation area at the top interface of current layer
184     real(r8)     :: conicw(pcols, kte)
185     real(r8)     :: cmfdqr(pcols, kte)
186     real(r8)     :: evapc(pcols, kte)         ! Evaporation rate of convective precipitation  
187     real(r8)     :: rainmr(pcols, kte)        ! rain mixing ratio
188     real(r8)     :: dlf(pcols,kte)            ! Detrainment of convective condensate
189     real(r8)     :: dlf2(pcols,kte)           ! Detrainment of convective condensate
190     real(r8)     :: cmfmc(pcols,kte+1)
191     real(r8)     :: cmfmc2(pcols,kte+1)
192     real(r8)     :: calday                    ! current calendar day
193     
194     
195     !State variables [The CAM's side variables are kept almost same by just 
196     !replacing '%' with an '_' (e.g state%t in CAM is state_t in WRF)]
197     real(r8)     :: state_t(pcols,kte)
198     real(r8)     :: state_q(pcols,kte,pcnst)
199     real(r8)     :: state_pmid(pcols,kte)
200     real(r8)     :: state_pdel(pcols,kte)
201     
202     real(r8)     :: ptend_s(pcols,kte)         !Dummy arguments for physics_update call
203     real(r8)     :: state_s(pcols,kte)         !Dummy arguments for physics_update call
204     
205     real(r8)     :: ptend_q(pcols,kte,pcnst)
206     
207     !Cam_out variable is a dummy placeholder as currently it is not used in parameterization
208     real(r8)     :: cam_out
210     
211     !! physics buffer 
212     ! Variables stored in physics buffer in CAM
213     real(r8), dimension(pcols,kte)       :: cldn       !cloud fraction
214     real(r8), dimension(pcols,kte)       :: cme        !local condensation of cloud water
215     real(r8), dimension(pcols,kte)       :: prain      !production of rain
216     real(r8), dimension(pcols,kte)       :: evapr      !evaporation of rain
217     real(r8), dimension(pcols,kte)       :: icwmrdp    !in cloud water mixing ratio, deep convection
218     real(r8), dimension(pcols,kte)       :: rprddp     !rain production, deep convection
219     real(r8), dimension(pcols,kte)       :: icwmrsh    !in cloud water mixing ratio, deep convection
220     real(r8), dimension(pcols,kte)       :: rprdsh     !rain production, deep convection
221     real(r8), dimension(pcols,kte,pcnst) :: fracis     !fraction of transported species that are insoluble
222     
223     !! Dec.29.2009. Sungsu
224     real(r8), dimension(pcols,kte)       ::  sh_frac   !Shallow convective cloud fraction
225     real(r8), dimension(pcols,kte)       ::  dp_frac   !Deep convective cloud fraction
226     real(r8), dimension(pcols,kte)       ::  evapcsh   !Evaporation rate of shallow convective precipitation >=0.
227     real(r8), dimension(pcols,kte)       ::  evapcdp   !Evaporation rate of deep    convective precipitation >=0.
228     !! Dec.29.2009. Sungsu 
229     
230 #if ( defined MODAL_AERO )
231     real(r8), dimension(pcols,kte,ntot_amode) :: dgnum_pbuf, dgnumwet_pbuf !Wet/ambient geom. mean diameter (m)
232     
233     !! for number distribution
234     real(r8), dimension(pcols,kte,pcnst)      :: qqcw     !cloud-borne aerosol
235     real(r8), dimension(pcols,kte,ntot_amode) :: qaerwat  !aerosol water
236     real(r8), dimension(pcols,kte)            :: rate1ord_cw2pr_st   !1st order rate for direct conversion of strat. cloud water to precip (1/s)    ! rce 2010/05/01
237 #endif
238     
239     nstep = itimestep
240     if(itimestep == 1) then
241        if(config_flags%shcu_physics .NE. CAMUWSHCUSCHEME) call wrf_message('WARNING: sh_frac,evapcsh,icwmrsh,rprdsh,cmfmc,cmfmc2  are set to zero in CAM_MAM_WETSCAV')
242        if(config_flags%cu_physics   .NE. CAMZMSCHEME)     call wrf_message('WARNING: evapcdp,icwmrdp,rprddp,dlf,dlf2 are set to zero in CAM_MAM_WETSCAV')
243     endif
244     !Initialize ptend_name,ptend_q,ptend_s,ptend_lq and ptend_ls so that ptend_q is zeroed out
245     call physics_ptend_init(ptend_name,ptend_q,ptend_s,ptend_lq,ptend_ls,pcnst)
246     
247     !Following arrays are declared as NaN as they are NOT used currently but still passed into the parameterization
248     cme(:,:) = nan
249     
250     !*NOTE*In CAM bcphiwet,ocphiwet,dstwet1,dstwet2,dstwet3 and dstwet4 of cam_out data struture are updated to be used for surface model.This has 
251     !NOT been implemented in WRF yet
252     cam_out = nan
253     
254     !Dry static energy is initialized to NaN and its tendency flag is set to .FALSE.
255     !Dry static enery is NOT required for wetscavaging but it is required to call WRF 
256     !implementation of CAM's physics_update in CAMMGMP microphysics
257     state_s(:,:) = nan
258     ptend_ls     = .FALSE.
259     ptend_s(:,:) = nan
260     
261     !Required assignments
262     p1st  = param_first_scalar ! Obtain CHEM array's first element's index
263     dt    = real(dtstep,r8)    ! Convert the time step to real*8
264     
265     ncol  = pcols
266     !This subroutine requires that ncol == 1
267     if(ncol .NE. 1) then
268        call wrf_error_fatal('Number of CAM Columns (NCOL) in CAM_MAM_WETSCAV scheme must be 1')
269     endif
270     
271     !Following varibales are set to NaN so that an error is produced whenever they are used
272     !calday is used only for Chemistry calculations inside the wetscav parameterization, which
273     !is not implemented yet in WRF
274     calday = nan
275     
276     !Divide domain in chuncks and map WRF variables into CAM
277     !Loop counters are named iw,jw,kw to represent that they index WRF sided variables
278     !The CAM's side variables are kept almost same by just replacing '%' with an '_' [e.g state%t in CAM is state_t in WRF]
279     
280     itsm1     = its - 1 
281     itile_len = ite - itsm1
282     do jw     = jts , jte 
283        do iw  = its , ite 
284           
285           lchnk   = (jw - jts) * itile_len + (iw - itsm1)             !1-D index location from a 2-D tile
286           ktep1   = kte + 1
287           
288           !Flip vertically quantities computed at the mid points        
289           do kw  = kts, kte
290              kflip                = ktep1 - kw
291              
292              state_pmid(1,kflip)  = p_hyd(iw,kw,jw)                   !Pressure at the mid-points (Pa) [state%pmid in CAM]  
293              dp                   = p8w(iw,kw,jw) - p8w(iw,kw+1,jw)   !Change in pressure (Pa) 
294              state_pdel(1,kflip)  = dp
295              state_t(1,kflip)     = t_phy(iw,kw,jw)                   !Temprature at the mid points (K) [state%t in CAM]
296              
297              !Following three formulas are obtained from ported CAM's ZM cumulus scheme
298              !Values of 0 cause a crash in entropy
299              multFrc              = 1._r8/(1._r8 + qv_curr(iw,kw,jw))
300              state_q(1,kflip,1)   = max( qv_curr(iw,kw,jw)*multFrc, 1.0e-30_r8 ) !Specific humidity                       [state%q(:,:,1) in CAM]
301              state_q(1,kflip,2)   = qc_curr(iw,kw,jw)*multFrc                    !Convert to moist mix ratio-cloud liquid [state%q(:,:,2) in CAM]
302              state_q(1,kflip,3)   = qi_curr(iw,kw,jw)*multFrc                    !cloud ice                               [state%q(:,:,3) in CAM]
303              state_q(1,kflip,4)   = nc3d(iw,kw,jw)*multFrc                       !Liquid cloud number
304              state_q(1,kflip,5)   = ni3d(iw,kw,jw)*multFrc                       !Ice cloud number
305              
306              !populate state_q and qqcw arrays
307              !Following Do-Loop is obtained from chem/module_cam_mam_aerchem_driver.F 
308              do l = p1st, num_chem
309                 l2 = lptr_chem_to_q(l)
310                 if ((l2 >= 1) .and. (l2 <= pcnst)) then
311                    state_q(1,kflip,l2) = chem(iw,kw,jw,l)*factconv_chem_to_q(l)
312                 end if
313                 l2 = lptr_chem_to_qqcw(l)
314                 if ((l2 >= 1) .and. (l2 <= pcnst)) then
315                    qqcw(1,kflip,l2) = chem(iw,kw,jw,l)*factconv_chem_to_q(l)     !Cloud borne aerosols
316                 end if
317              end do ! l
318              
319              !Populate dgnums appropriately
320              !Following Do-Loop is obtained from chem/module_cam_mam_aerchem_driver.F 
321              do imode = 1 , ntot_amode
322                 dgnum_pbuf(1,kflip,imode)    = dgnum4D(iw,kw,jw,imode)           !Obtained from 4D arrays 
323                 dgnumwet_pbuf(1,kflip,imode) = dgnumwet4D(iw,kw,jw,imode) 
324                 
325                 l = waterptr_aer(1,imode)
326                 if ((l >= p1st) .and. (l <= num_chem)) then
327                    qaerwat(1,kflip,imode) = chem(iw,kw,jw,l)*factconv_chem_to_q(l)!aerosol water 
328                 endif
329              enddo
330              
331              !*NOTE* QME3D doesn't include contribution from MACROPHYSICS (CMELIQ).Assinment to CME is
332              !Commented out currently as CME is NEVER used in wetscaging code
333              !cme(1,kflip)         = qme3d(iw,kw,jw)
334              cme(1,kflip)         = nan                                          !Net condensation rate (kg/kg/s)
335              prain(1,kflip)       = prain3d(iw,kw,jw)                            !Rate of conversion of condensate to precipitation (kg/kg/s)
336              evapr(1,kflip)       = nevapr3d(iw,kw,jw)                           !Evaporation rate of rain + snow (kg/kg/s)
337              
338              rate1ord_cw2pr_st(1,kflip) = rate1ord_cw2pr_st3d(iw,kw,jw)          !1st order rate for direct conversion of strat. cloud water to precip (1/s)
339              if(is_CAMMGMP_used) then
340                 cldn(1,kflip)              = cldfra_mp_all(iw,kw,jw)                !Cloud fraction
341              else
342                 cldn(1,kflip)              = cldfra(iw,kw,jw)
343              endif
344              cldn(1,kflip)              = min(max(cldn(1,kflip),0._r8),1._r8)
345              
346              sh_frac(1,kflip)           = 0.0_r8
347              evapcsh(1,kflip)           = 0.0_r8
348              icwmrsh(1,kflip)           = 0.0_r8
349              rprdsh(1,kflip)            = 0.0_r8
350              
351              if(config_flags%shcu_physics==CAMUWSHCUSCHEME) then
352                 !inputs from shallow convection
353                 sh_frac(1,kflip)        = shfrc3d(iw,kw,jw)                      !Shallow cloud fraction         
354                 evapcsh(1,kflip)        = evapcsh3d(iw,kw,jw)                    !Evaporation of shallow convection precipitation (kg/kg/s)
355                 icwmrsh(1,kflip)        = icwmrsh3d(iw,kw,jw)                    !shallow cumulus in-cloud water mixing ratio (kg/m2)
356                 rprdsh(1,kflip)         = rprdsh3d(iw,kw,jw)                     !dq/dt due to deep and shallow convective rainout(kg/kg/s)
357              endif
358              
359              evapcdp(1,kflip)           = 0.0_r8
360              icwmrdp(1,kflip)           = 0.0_r8
361              rprddp(1,kflip)            = 0.0_r8
362              dlf(1,kflip)               = 0.0_r8
363              dlf2(1,kflip)              = 0.0_r8
364              
365              if(config_flags%cu_physics==CAMZMSCHEME)then
366                 !inputs from deep convection
367                 evapcdp(1,kflip)        = evapcdp3d(iw,kw,jw)                    !Evaporation of deep convective precipitation (kg/kg/s)
368                 icwmrdp(1,kflip)        = icwmrdp3d(iw,kw,jw)                    !Deep Convection in-cloud water mixing ratio (kg/m2)
369                 rprddp(1,kflip)         = rprddp3d(iw,kw,jw)                     !dq/dt due to deep convective rainout (kg/kg/s)
370                 dlf(1,kflip)            = dlf3d(iw,kw,jw)                        !Detrainment of convective condensate (kg/kg/s)
371                 dlf2(1,kflip)           = dlf2_3d(iw,kw,jw)                      !dq/dt due to export of cloud water into environment by shallow convection(kg/kg/s)  
372              endif
373           enddo
374           
375           do kw = kts, kte+1
376              kflip = kte - kw + 2
377              
378              cmfmc(1,kflip)      = 0.0_r8
379              cmfmc2(1,kflip)     = 0.0_r8
380              if(config_flags%shcu_physics==CAMUWSHCUSCHEME) then
381                 cmfmc(1,kflip)   = cmfmc3d(iw,kw,jw)    !Deep + Shallow Convective mass flux [ kg /s/m^2 ]
382                 cmfmc2(1,kflip)  = cmfmc2_3d(iw,kw,jw)  !Shallow convective mass flux [ kg/s/m^2 ]
383              endif
384           end do
385           
386           do kcam = 1, kte
387              !Formulation for dp_frac is obtained from cloud_fraction.F90 of CAM
388              dp_frac(1,kcam)         = max(0.0_r8,min(dp1*log(1.0_r8+dp2*(cmfmc(1,kcam+1)-cmfmc2(1,kcam+1))),0.60_r8))
389           end do
390           
391           !The CAM wet scavenging computations begin here
392           cldc(:ncol,:)  = dp_frac(:ncol,:) + sh_frac(:ncol,:) !! Sungsu included this.
393           evapc(:ncol,:) = evapcsh(:ncol,:) + evapcdp(:ncol,:) !! Sungsu included this.
394           clds(:ncol,:)  = cldn(:ncol,:) - cldc(:ncol,:)       !! Stratiform cloud fraction
395           
396           
397           !! sum deep and shallow convection contributions
398           conicw(:ncol,:) = (icwmrdp(:ncol,:)*dp_frac(:ncol,:) + icwmrsh(:ncol,:)*sh_frac(:ncol,:))/ &
399                max(0.01_r8, sh_frac(:ncol,:) + dp_frac(:ncol,:))
400           
401           cmfdqr(:ncol,:) = rprddp(:ncol,:)  + rprdsh(:ncol,:)
402           
403           
404           !OUTPUT- cldv, cldvcu, cldvst and rainmr 
405           
406           !!   fields needed for wet scavenging
407           call clddiag( state_t, state_pmid, state_pdel, cmfdqr, evapc, cldn, cldc, clds, cme, evapr, prain, &
408                cldv, cldvcu, cldvst, rainmr, ncol )
409           
410           ptend_name = 'wetdep'
411           
412           !*Please Note:* Calls to modal_aero_calcsize_sub and modal_aero_wateruptake_sub are taken care of in module_cam_mam_aerchem_driver.F
413           
414           !Output- ptend_name,ptend_lq,ptend_q, fracis, qqcw, qaerwat
415           
416           !Balwinder.Singh@pnnl.gov: Changed the arguments to the following 
417           ! call in CAM so that 'state','ptend' and 'cam_out' data structures are not 
418           ! passed into the call.
419           fracis(:,:,:) = 1.0_r8
420           call mz_aero_wet_intr (lchnk, ncol, state_q,                 &
421                state_pdel, state_pmid, state_t, ptend_name,            &
422                ptend_lq, ptend_q, nstep, dt, cme, prain, evapr, cldv,  &
423                cldvcu, cldvst, cldc, cldn, fracis, calday, cmfdqr,     &
424                evapc, conicw, rainmr,                                  &
425                rate1ord_cw2pr_st,                                      &   ! rce 2010/05/01
426                dgnumwet_pbuf, qqcw, qaerwat, cam_out, dlf              )
428           call physics_update(lchnk,dt,state_q,ptend_q,state_s,ptend_s,ptend_name,ptend_lq,ptend_ls, pcnst)
429           
430           !Post processing of the output from CAM's parameterization
431           do kw=kts,kte
432              kflip = kte-kw+1
433              do imode = 1,  ntot_amode
434                 l = waterptr_aer(1,imode)
435                 if ((l >= p1st) .and. (l <= num_chem)) then
436                    chem(iw,kw,jw,l) = qaerwat(1,kflip,imode)/factconv_chem_to_q(l)
437                 endif
438              end do ! imode
439              
440              !Following equation are derived following UWPBL and CAMZM schemes
441              qv_curr(iw,kw,jw)       = state_q(1,kflip,1) / (1.0_r8 - state_q(1,kflip,1)) 
442              multFrc                 = 1._r8 + qv_curr(iw,kw,jw)
443              
444              qc_curr(iw,kw,jw)       = state_q(1,kflip,2) * multFrc
445              qi_curr(iw,kw,jw)       = state_q(1,kflip,3) * multFrc 
446              nc3d(iw,kw,jw)          = state_q(1,kflip,4) * multFrc  
447              ni3d(iw,kw,jw)          = state_q(1,kflip,5) * multFrc
448              do l = 1 ,5
449                 fracis3d(iw,kw,jw,l)     = fracis(1,kflip,l)          !Fraction of transported species that are insoluble             
450              enddo
451              do l = p1st, num_chem
452                 l2 = lptr_chem_to_q(l)
453                 if ((l2 >= 1) .and. (l2 <= pcnst)) then
454                    chem(iw,kw,jw,l) = state_q(1,kflip,l2)/factconv_chem_to_q(l)
455                    fracis3d(iw,kw,jw,l2)     = fracis(1,kflip,l2)          !Fraction of transported species that are insoluble             
456                 end if
457                 l2 = lptr_chem_to_qqcw(l)
458                 if ((l2 >= 1) .and. (l2 <= pcnst)) then
459                    chem(iw,kw,jw,l) = qqcw(1,kflip,l2)/factconv_chem_to_q(l)
460                 end if
461              end do ! l
462           end do
463           
464        enddo !iw loop
465     enddo !jw loop
466     return
467     
468   end subroutine wetscav_cam_mam_driver
469   
470   !----------------------------------------------------------------------------
471   subroutine wetscav_cam_mam_driver_init(ids,ide, jds,jde, kds,kde, &
472        ims,ime, jms,jme, kms,kme,                           &
473        its,ite, jts,jte, kts,kte                            )
474     !    
475     !Purpose: 
476     !Initialize few variables needed for the driver
477     !     
478     !Author: Balwinder.Singh@pnnl.gov
479     !----------------------------------------------------------------------------
480     use module_cam_support,        only: pver
481     use mz_aerosols_intr,          only: modal_aero_bcscavcoef_init, mz_aero_initialize 
482     implicit none
483     integer, intent(in) :: ids,ide, jds,jde, kds,kde
484     integer, intent(in) :: ims,ime, jms,jme, kms,kme
485     integer, intent(in) :: its,ite, jts,jte, kts,kte
486     
487     jtf   = min(jte,jde-1)
488     ktf   = min(kte,kde-1)
489     itf   = min(ite,ide-1)
491     !Map CAM veritcal level variables
492     pver  = ktf - kts + 1 
493     pverp = pver + 1
495     !Following constants (dp1 and dp2) are  obtained from cloud_fraction.F90 of CAM for highest resolution(0.23x0.31)
496     dp1   = 0.10_r8 
497     dp2   = 500.0_r8
498     
499     !initialize mz_aero
500     call  mz_aero_initialize 
501     
502   end subroutine wetscav_cam_mam_driver_init
503   
504 end module module_cam_mam_wetscav