1 MODULE module_dep_simple
4 ! many of these parameters will depend on the RADM mechanism!
5 ! if you change it, lets talk about it and get it done!!!
8 INTEGER, PARAMETER :: dep_seasons = 5 , nlu = 25, &
9 nseason = 1, nseasons = 2
11 ! following currently hardwired to USGS
13 INTEGER, PARAMETER :: isice_temp=24,iswater_temp=16
14 character, parameter :: mminlu='USGS'
16 INTEGER :: ixxxlu(nlu)
18 rac(nlu,dep_seasons), rclo(nlu,dep_seasons), rcls(nlu,dep_seasons), &
19 rgso(nlu,dep_seasons), rgss(nlu,dep_seasons), &
20 ri(nlu,dep_seasons), rlu(nlu,dep_seasons)
22 ! NO MORE THAN 1000 SPECIES FOR DEPOSITION
24 REAL, DIMENSION (1:1000) :: dratio,hstar,hstar4,f0,dhr,scpr23
25 ! .. Default Accessibility ..
27 logical, allocatable :: is_aerosol(:) ! true if field is aerosol (any phase)
30 subroutine wesely_driver(id,ktau,dtstep, &
32 gmt,julday,t_phy,moist,p8w,t8w,raincv, &
33 p_phy,chem,rho_phy,dz8w,ddvel,aer_res_def,aer_res_zcen, &
34 ivgtyp,tsk,gsw,vegfra,pbl,rmol,ust,znt,xlat,xlong,z,z_at_w,&
36 ids,ide, jds,jde, kds,kde, &
37 ims,ime, jms,jme, kms,kme, &
38 its,ite, jts,jte, kts,kte )
39 !----------------------------------------------------------------------
40 USE module_model_constants
42 USE module_state_description
43 USE module_data_sorgam
45 INTEGER, INTENT(IN ) :: id,julday, &
47 ids,ide, jds,jde, kds,kde, &
48 ims,ime, jms,jme, kms,kme, &
49 its,ite, jts,jte, kts,kte
50 INTEGER, INTENT(IN ) :: &
52 REAL, INTENT(IN ) :: &
55 ! advected moisture variables
57 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_moist ), &
60 ! advected chemical species
62 REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ), &
63 INTENT(INOUT ) :: chem
65 ! deposition velocities
67 REAL, DIMENSION( its:ite, jts:jte, num_chem ), &
68 INTENT(INOUT ) :: ddvel
70 ! input from met model
73 REAL, DIMENSION( ims:ime , kms:kme , jms:jme ) , &
81 INTEGER,DIMENSION( ims:ime , jms:jme ) , &
84 REAL, DIMENSION( ims:ime , jms:jme ) , &
96 !--- deposition and emissions stuff
99 REAL :: clwchem,dvfog,dvpart,pa,rad,rhchem,ta,ustar,vegfrac,z1,zntt
100 INTEGER :: i, iland, iprt, iseason, j, jce, jcs, n, nr, ipr,jpr,nvr
101 LOGICAL :: highnh3, rainflag, vegflag, wetflag
104 REAL :: p(kts:kte), srfres(numgas),ddvel0d(numgas)
106 ! necessary for aerosols (module dependent)
108 real :: aer_res_def(its:ite,jts:jte), aer_res_zcen(its:ite,jts:jte)
111 TYPE(grid_config_rec_type), INTENT(IN ) :: config_flags
115 ! .. Intrinsic Functions ..
118 CALL wrf_debug(15,'in dry_dep_wesely')
120 if(julday.lt.90.or.julday.gt.270)then
122 CALL wrf_debug(15,'setting iseason to 2')
130 vegfrac = vegfra(i,j)
131 pa = .01*p_phy(i,kts,j)
132 clwchem = moist(i,kts,j,p_qc)
136 z1 = z_at_w(i,kts+1,j)-z_at_w(i,kts,j)
138 ! Set logical default values
143 if(moist(i,kts,j,p_qr).gt.1.e-18 .or. raincv(i,j).gt.0.)rainflag = .true.
145 rhchem = MIN( 100.,100. * moist(i,kts,j,p_qv) / &
146 (3.80*exp(17.27*(t_phy(i,kts,j)-273.)/(t_phy(i,kts,j)-36.))/pa))
147 rhchem = MAX(5.,RHCHEM)
148 if (rhchem >= 95.) wetflag = .true.
151 if(chem(i,kts,j,p_nh3).gt.2.*chem(i,kts,j,p_so2))highnh3 = .true.
156 ! if(snowc(i,j).gt.0.)iseason=4
157 CALL rc(rcx,ta,rad,rhchem,iland,iseason,numgas, &
158 wetflag,rainflag,highnh3,iprt,p_o3,p_so2,p_nh3)
163 CALL deppart(rmol(i,j),ustar,rhchem,clwchem,iland,dvpart,dvfog)
167 CALL landusevg(ddvel0d,ustar,rmol(i,j),zntt,z1,dvpart,iland, &
168 numgas,srfres,aer_res_def(i,j),aer_res_zcen(i,j),p_sulf)
170 !wig: CMBZ does not have HO and HO2 last so need to copy all species
171 ! ddvel(i,j,1:numgas-2)=ddvel0d(1:numgas-2)
172 ddvel(i,j,1:numgas)=ddvel0d(1:numgas)
173 if ( (config_flags%chem_opt == RADM2 ) .or. &
174 (config_flags%chem_opt == RADM2SORG ) .or. &
175 (config_flags%chem_opt == RADM2SORG_AQ) ) then
176 ddvel(i,j,p_hcl) = ddvel(i,j,p_hno3)
182 ! For the additional CBMZ species, assign similar RADM counter parts for
183 ! now. Short lived species get a zero velocity since dry dep should be
184 ! unimportant. **ALSO**, treat p_sulf as h2so4 vapor, not aerosol sulfate
186 if ( (config_flags%chem_opt == CBMZ ) .or. &
187 (config_flags%chem_opt == CBMZ_BB ) .or. &
188 (config_flags%chem_opt == CBMZ_BB_KPP ) .or. &
189 (config_flags%chem_opt == CBMZ_MOSAIC_4BIN_AQ) .or. &
190 (config_flags%chem_opt == CBMZ_MOSAIC_8BIN_AQ) .or. &
191 (config_flags%chem_opt == CBMZ_MOSAIC_4BIN) .or. &
192 (config_flags%chem_opt == CBMZ_MOSAIC_8BIN) ) then
195 ddvel(i,j,p_sulf) = ddvel(i,j,p_hno3)
196 ddvel(i,j,p_hcl) = ddvel(i,j,p_hno3)
197 ddvel(i,j,p_ch3o2) = 0
198 ddvel(i,j,p_ethp) = 0
199 ddvel(i,j,p_ch3oh) = ddvel(i,j,p_hcho)
200 ddvel(i,j,p_c2h5oh) = ddvel(i,j,p_hcho)
201 !wig, 1-May-2007 (added par to tables) ddvel(i,j,p_par) = ddvel(i,j,p_hc5)
204 ddvel(i,j,p_open) = ddvel(i,j,p_xyl)
205 ddvel(i,j,p_op3) = ddvel(i,j,p_op2)
206 ddvel(i,j,p_c2o3) = 0
208 ddvel(i,j,p_ano2) = 0
211 ddvel(i,j,p_xpar) = 0
212 ddvel(i,j,p_isoprd) = 0
213 ddvel(i,j,p_isopp) = 0
214 ddvel(i,j,p_isopn) = 0
215 ddvel(i,j,p_isopo2) = 0
216 if( config_flags%chem_opt == CBMZ ) then
218 ddvel(i,j,p_msa) = ddvel(i,j,p_hno3)
219 ddvel(i,j,p_dmso) = 0
220 ddvel(i,j,p_dmso2) = 0
221 ddvel(i,j,p_ch3so2h) = 0
222 ddvel(i,j,p_ch3sch2oo) = 0
223 ddvel(i,j,p_ch3so2) = 0
224 ddvel(i,j,p_ch3so3) = 0
225 ddvel(i,j,p_ch3so2oo) = 0
226 ddvel(i,j,p_ch3so2ch2oo) = 0
233 ! For the additional CBM4 species, assign similar RADM counter parts for
234 ! now. Short lived species get a zero velocity since dry dep should be
237 if (config_flags%chem_opt == CBM4_KPP ) then
240 ddvel(i,j,p_open) = ddvel(i,j,p_xyl)
243 ddvel(i,j,p_ald2) = ddvel(i,j,p_ald)
249 ! For gocartracm,radm
251 if ((config_flags%chem_opt == GOCARTRACM_KPP) .OR. &
252 (config_flags%chem_opt == GOCARTRADM2_KPP) .OR. &
253 (config_flags%chem_opt == GOCARTRADM2)) then
256 ddvel(i,j,p_sulf) = 0.
257 ddvel(i,j,p_dms) = 0.
258 ddvel(i,j,p_msa) = ddvel(i,j,p_hno3)
259 if( config_flags%chem_opt == GOCARTRADM2 ) then
260 ddvel(i,j,p_hcl) = ddvel(i,j,p_hno3)
265 ! For gocartsimple : need msa. On the other hand sulf comes from aerosol routine
267 if (config_flags%chem_opt == GOCART_SIMPLE ) then
270 ddvel(i,j,p_msa) = ddvel(i,j,p_sulf)
271 ddvel(i,j,p_sulf) = 0.
272 ddvel(i,j,p_dms) = 0.
277 END SUBROUTINE wesely_driver
279 ! **********************************************************************
280 ! ************************** SUBROUTINE RC ***************************
281 ! **********************************************************************
282 SUBROUTINE rc(rcx,t,rad,rh,iland,iseason,numgas, &
283 wetflag,rainflag,highnh3,iprt,p_o3,p_so2,p_nh3)
284 ! THIS SUBROUTINE CALCULATES SURFACE RESISTENCES ACCORDING
287 ! ATMOSPHERIC ENVIRONMENT 23 (1989), 1293-1304
288 ! WITH SOME ADDITIONS ACCORDING TO
289 ! J. W. ERISMAN, A. VAN PUL, AND P. WYERS,
290 ! ATMOSPHERIC ENVIRONMENT 28 (1994), 2595-2607
291 ! WRITTEN BY WINFRIED SEIDL, APRIL 1997
292 ! MODYFIED BY WINFRIED SEIDL, MARCH 2000
294 !----------------------------------------------------------------------
295 ! .. Scalar Arguments ..
297 INTEGER :: iland, iseason, numgas
298 LOGICAL :: highnh3, rainflag, wetflag
301 INTEGER :: iprt,p_o3,p_so2,p_nh3
303 ! .. Local Scalars ..
304 REAL :: rclx, rdc, resice, rgsx, rluo1, rluo2, rlux, rmx, rs, rsmx, &
309 REAL :: hstary(numgas)
311 ! .. Intrinsic Functions ..
323 !!! HARDWIRE VALUES FOR TESTING
331 IF ((tc<=0.) .OR. (tc>=40.)) THEN
334 rs = ri(iland,iseason)*(1+z*z)*(400./(tc*(40.-tc)))
336 rdc = 100*(1.+1000./(rad+10))/(1+1000.*rdtheta)
337 rluo1 = 1./(1./3000.+1./3./rlu(iland,iseason))
338 rluo2 = 1./(1./1000.+1./3./rlu(iland,iseason))
339 resice = 1000.*exp(-tc-4.)
342 IF (hstar(n)==0.) GO TO 10
343 hstary(n) = hstar(n)*exp(dhr(n)*(1./t-1./298.))
344 rmx = 1./(hstary(n)/3000.+100.*f0(n))
345 rsmx = rs*dratio(n) + rmx
346 rclx = 1./(hstary(n)/1.E+5/rcls(iland,iseason)+f0(n)/rclo(iland, &
348 rgsx = 1./(hstary(n)/1.E+5/rgss(iland,iseason)+f0(n)/rgso(iland, &
350 rlux = rlu(iland,iseason)/(1.E-5*hstary(n)+f0(n)) + resice
352 rlux = 1./(1./3./rlu(iland,iseason)+1.E-7*hstary(n)+f0(n)/rluo1)
355 rlux = 1./(1./3./rlu(iland,iseason)+1.E-7*hstary(n)+f0(n)/rluo2)
357 rcx(n) = 1./(1./rsmx+1./rlux+1./(rdc+rclx)+1./(rac(iland, &
359 IF (rcx(n)<1.) rcx(n) = 1.
362 ! SPECIAL TREATMENT FOR OZONE
364 hstary(p_o3) = hstar(p_o3)*exp(dhr(p_o3)*(1./t-1./298.))
365 rmx = 1./(hstary(p_o3)/3000.+100.*f0(p_o3))
366 rsmx = rs*dratio(p_o3) + rmx
367 rlux = rlu(iland,iseason)/(1.E-5*hstary(p_o3)+f0(p_o3)) + resice
368 rclx = rclo(iland,iseason) + resice
369 rgsx = rgso(iland,iseason) + resice
370 IF (wetflag) rlux = rluo1
371 IF (rainflag) rlux = rluo2
372 rcx(p_o3) = 1./(1./rsmx+1./rlux+1./(rdc+rclx)+1./(min(100.,rac(iland, &
374 IF (rcx(p_o3)<1.) rcx(p_o3) = 1.
377 ! SPECIAL TREATMENT FOR SO2 (Wesely)
378 ! HSTARY(P_SO2)=HSTAR(P_SO2)*EXP(DHR(P_SO2)*(1./T-1./298.))
379 ! RMX=1./(HSTARY(P_SO2)/3000.+100.*F0(P_SO2))
380 ! RSMX=RS*DRATIO(P_SO2)+RMX
381 ! RLUX=RLU(ILAND,ISEASON)/(1.E-5*HSTARY(P_SO2)+F0(P_SO2))
383 ! RCLX=RCLS(ILAND,ISEASON)+RESICE
384 ! RGSX=RGSS(ILAND,ISEASON)+RESICE
385 ! IF ((wetflag).OR.(RAINFLAG)) THEN
386 ! IF (ILAND.EQ.1) THEN
392 ! RCX(P_SO2)=1./(1./RSMX+1./RLUX+1./(RDC+RCLX)
393 ! & +1./(RAC(ILAND,ISEASON)+RGSX))
394 ! IF (RCX(P_SO2).LT.1.) RCX(P_SO2)=1.
396 ! SO2 according to Erisman et al. 1994
398 rsmx = rs*dratio(p_so2)
402 rlux = 25000.*exp(-0.0693*rh)
404 rlux = 0.58E12*exp(-0.278*rh)
407 IF (((wetflag) .OR. (rainflag)) .AND. (tc>(-1.))) THEN
410 IF ((tc>=(-5.)) .AND. (tc<=(-1.))) THEN
416 ! INSTEAD OF R_INC R_CL and R_DC of Wesely are used
417 rclx = rcls(iland,iseason)
421 IF ((wetflag) .OR. (rainflag)) THEN
429 IF (iland==iswater_temp) THEN
433 IF ((iseason==4) .OR. (iland==isice_temp)) THEN
437 IF ((tc>=(-1.)) .AND. (tc<=2.)) THEN
444 ! TOTAL SURFACE RESISTENCE
445 IF ((iseason/=4) .AND. (ixxxlu(iland)/=1) .AND. (iland/=iswater_temp) .AND. &
446 (iland/=isice_temp)) THEN
447 rcx(p_so2) = 1./(1./rsmx+1./rlux+1./(rclx+rdc+rgsx))
451 IF (rcx(p_so2)<1.) rcx(p_so2) = 1.
452 ! NH3 according to Erisman et al. 1994
459 rsmx = rs*dratio(p_nh3)
460 ! GRASSLAND (PASTURE DURING GRAZING)
461 IF (ixxxlu(iland)==3) THEN
466 IF ((iseason==2) .OR. (iseason==3) .OR. (iseason==5)) THEN
474 IF ((wetflag) .OR. (rainflag)) THEN
478 IF ((tc>=(-5.)) .AND. (tc<=-1.)) THEN
486 ! AGRICULTURAL LAND (CROPS AND UNGRAZED PASTURE)
487 IF (ixxxlu(iland)==2) THEN
495 IF ((wetflag) .OR. (rainflag)) THEN
499 IF ((iseason==2) .OR. (iseason==3) .OR. (iseason==5)) THEN
507 IF ((wetflag) .OR. (rainflag)) THEN
511 IF ((tc>=(-5.)) .AND. (tc<=-1.)) THEN
519 ! SEMI-NATURAL ECOSYSTEMS AND FORESTS
520 IF ((ixxxlu(iland)==4) .OR. (ixxxlu(iland)==5) .OR. (ixxxlu( &
527 IF ((wetflag) .OR. (rainflag)) THEN
534 IF ((iseason==2) .OR. (iseason==3) .OR. (iseason==5)) THEN
536 IF ((tc>=(-5.)) .AND. (tc<=-1.)) THEN
545 IF (iland==iswater_temp) THEN
548 ! URBAN AND DESERT (SOIL SURFACES)
549 IF (ixxxlu(iland)==1) THEN
550 IF ( .NOT. wetflag) THEN
556 ! SNOW COVERED SURFACES OR PERMANENT ICE
557 IF ((iseason==4) .OR. (iland==isice_temp)) THEN
561 IF ((tc>=(-1.)) .AND. (tc<=2.)) THEN
562 rcx(p_nh3) = 70.*(2.-tc)
568 IF (rcx(p_nh3)<1.) rcx(p_nh3) = 1.
571 ! **********************************************************************
572 ! ************************ SUBROUTINE DEPPART **************************
573 ! **********************************************************************
574 SUBROUTINE deppart(rmol,ustar,rh,clw,iland,dvpart,dvfog)
575 ! THIS SUBROUTINE CALCULATES SURFACE DEPOSITION VELOCITIES
576 ! FOR FINE AEROSOL PARTICLES ACCORDING TO THE MODEL OF
577 ! J. W. ERISMAN, A. VAN PUL, AND P. WYERS,
578 ! ATMOSPHERIC ENVIRONMENT 28 (1994), 2595-2607
579 ! WRITTEN BY WINFRIED SEIDL, APRIL 1997
580 ! MODIFIED BY WINFRIED SEIDL, MARCH 2000
582 ! ----------------------------------------------------------------------
583 ! .. Scalar Arguments ..
584 REAL :: clw, dvfog, dvpart, rh, rmol, ustar
587 ! .. Intrinsic Functions ..
590 dvpart = ustar/kpart(iland)
592 ! INSTABLE LAYERING CORRECTION
593 dvpart = dvpart*(1.+(-300.*rmol)**0.66667)
596 ! HIGH RELATIVE HUMIDITY CORRECTION
597 ! ACCORDING TO J. W. ERISMAN ET AL.
598 ! ATMOSPHERIC ENVIRONMENT 31 (1997), 321-332
599 dvpart = dvpart*(1.+0.37*exp((rh-80.)/20.))
602 ! SEDIMENTATION VELOCITY OF FOG WATER ACCORDING TO
603 ! R. FORKEL, W. SEIDL, R. DLUGI AND E. DEIGELE
604 ! J. GEOPHYS. RES. 95D (1990), 18501-18515
606 IF (ixxxlu(iland)==5) THEN
607 ! TURBULENT DEPOSITION OF FOG WATER IN CONIFEROUS FOREST ACCORDI
608 ! A. T. VERMEULEN ET AL.
609 ! ATMOSPHERIC ENVIRONMENT 31 (1997), 375-386
610 dvfog = dvfog + 0.195*ustar*ustar
613 END SUBROUTINE deppart
614 SUBROUTINE landusevg(vgs,ustar,rmol,z0,zz,dvparx,iland,numgas, &
615 srfres,aer_res_def,aer_res_zcen,p_sulf)
616 ! This subroutine calculates the species specific deposition velocit
617 ! as a function of the local meteorology and land use. The depositi
618 ! Velocity is also landuse specific.
619 ! Reference: Hsieh, C.M., Wesely, M.L. and Walcek, C.J. (1986)
620 ! A Dry Deposition Module for Regional Acid Deposition
621 ! EPA report under agreement DW89930060-01
622 ! Revised version by Darrell Winner (January 1991)
623 ! Environmental Engineering Science 138-78
624 ! California Institute of Technology
626 ! Modified by Winfried Seidl (August 1997)
627 ! Fraunhofer-Institut fuer Atmosphaerische Umweltforschung
628 ! Garmisch-Partenkirchen, D-82467
629 ! for use of Wesely and Erisman surface resistances
631 ! Ustar : The grid average friction velocity (m/s)
632 ! Rmol : Reciprocal of the Monin-Obukhov length (1/m)
633 ! Z0 : Surface roughness height for the grid square (m)
634 ! SrfRes : Array of landuse/atmospheric/species resistances (s/m)
635 ! Slist : Array of chemical species codes
636 ! Dvparx : Array of surface deposition velocity of fine aerosol p
638 ! Vgs : Array of species and landuse specific deposition
640 ! Vg : Cell-average deposition velocity by species (m/s)
642 ! SCPR23 : (Schmidt #/Prandtl #)**(2/3) Diffusion correction fac
643 ! Zr : Reference Height (m)
644 ! Iatmo : Parameter specifying the stabilty class (Function of
645 ! Z0 : Surface roughness height (m)
646 ! karman : Von Karman constant (from module_model_constants)
647 USE module_model_constants, only: karman
649 ! .. Scalar Arguments ..
650 REAL :: dvparx, rmol, ustar, z0, zz
651 real :: aer_res_def, aer_res_zcen, polint
652 INTEGER :: iland, numgas, p_sulf
654 ! .. Array Arguments ..
655 REAL :: srfres(numgas), vgs(numgas)
657 ! .. Local Scalars ..
658 REAL :: rmol_tmp, vgp, vgpart, zr
662 REAL :: vgspec(numgas)
664 ! Calculate aerodynamic resistance for reference height = layer center
667 CALL depvel(numgas,rmol_tmp,zr,z0,ustar,vgspec,vgpart,aer_res_zcen)
668 ! Set the reference height (10.0 m)
672 ! CALCULATE THE DEPOSITION VELOCITY without any surface
673 ! resistance term, i.e. 1 / (ra + rb)
674 CALL depvel(numgas,rmol,zr,z0,ustar,vgspec,vgpart,aer_res_def)
676 ! Calculate the deposition velocity for each species
677 ! and grid cell by looping through all the possibile combinations
680 vgp = 1.0/((1.0/vgpart)+(1.0/dvparx))
682 ! Loop through the various species
686 ! Add in the surface resistance term, rc (SrfRes)
688 vgs(jspec) = 1.0/(1.0/vgspec(jspec)+srfres(jspec))
692 CALL cellvg(vgs,ustar,zz,zr,rmol,numgas)
695 END SUBROUTINE landusevg
697 SUBROUTINE cellvg(vgtemp,ustar,dz,zr,rmol,nspec)
698 ! THIS PROGRAM HAS BEEN DESIGNED TO CALCULATE THE CELL AVERAGE
699 ! DEPOSITION VELOCITY GIVEN THE VALUE OF VG AT SOME REFERENCE
700 ! HEIGHT ZR WHICH IS MUCH SMALLER THAN THE CELL HEIGHT DZ.
701 ! PROGRAM WRITTEN BY GREGORY J.MCRAE (NOVEMBER 1977)
702 ! Modified by Darrell A. Winner (February 1991)
703 !.....PROGRAM VARIABLES...
704 ! VgTemp - DEPOSITION VELOCITY AT THE REFERENCE HEIGHT
705 ! USTAR - FRICTION VELOCITY
706 ! RMOL - RECIPROCAL OF THE MONIN-OBUKHOV LENGTH
707 ! ZR - REFERENCE HEIGHT
709 ! CELLVG - CELL AVERAGE DEPOSITION VELOCITY
710 ! VK - VON KARMAN CONSTANT
711 USE module_model_constants, only: karman
713 ! .. Scalar Arguments ..
714 REAL :: dz, rmol, ustar, zr
717 ! .. Array Arguments ..
718 REAL :: vgtemp(nspec)
720 ! .. Local Scalars ..
721 REAL :: a, fac, pdz, pzr, vk
724 ! .. Intrinsic Functions ..
727 ! Set the von Karman constant
730 ! DETERMINE THE STABILITY BASED ON THE CONDITIONS
738 pdz = sqrt(1.0-9.0*dz*rmol)
739 pzr = sqrt(1.0-9.0*zr*rmol)
740 fac = ((pdz-1.0)/(pzr-1.0))*((pzr+1.0)/(pdz+1.0))
741 a = 0.74*dz*alog(fac) + (0.164/rmol)*(pdz-pzr)
742 ELSE IF (rmol==0) THEN
743 a = 0.74*(dz*alog(dz/zr)-dz+zr)
745 a = 0.74*(dz*alog(dz/zr)-dz+zr) + (2.35*rmol)*(dz-zr)**2
748 ! CALCULATE THE DEPOSITION VELOCITIY
750 vgtemp(nss) = vgtemp(nss)/(1.0+vgtemp(nss)*a/(vk*ustar*(dz-zr)))
754 END SUBROUTINE cellvg
755 SUBROUTINE depvel(numgas,rmol,zr,z0,ustar,depv,vgpart,aer_res)
756 ! THIS FUNCTION HAS BEEN DESIGNED TO EVALUATE AN UPPER LIMIT
757 ! FOR THE POLLUTANT DEPOSITION VELOCITY AS A FUNCTION OF THE
758 ! SURFACE ROUGHNESS AND METEOROLOGICAL CONDITIONS.
759 ! PROGRAM WRITTEN BY GREGORY J.MCRAE (NOVEMBER 1977)
760 ! Modified by Darrell A. Winner (Feb. 1991)
761 ! by Winfried Seidl (Aug. 1997)
762 !.....PROGRAM VARIABLES...
763 ! RMOL - RECIPROCAL OF THE MONIN-OBUKHOV LENGTH
764 ! ZR - REFERENCE HEIGHT
765 ! Z0 - SURFACE ROUGHNESS HEIGHT
766 ! SCPR23 - (Schmidt #/Prandtl #)**(2/3) Diffusion correction fact
767 ! UBAR - ABSOLUTE VALUE OF SURFACE WIND SPEED
768 ! DEPVEL - POLLUTANT DEPOSITION VELOCITY
769 ! Vk - VON KARMAN CONSTANT
770 ! USTAR - FRICTION VELOCITY U*
771 ! POLINT - POLLUTANT INTEGRAL
772 ! AER_RES - AERODYNAMIC RESISTANCE
774 ! MCRAE, G.J. ET AL. (1983) MATHEMATICAL MODELING OF PHOTOCHEMICAL
775 ! AIR POLLUTION, ENVIRONMENTAL QUALITY LABORATORY REPORT 18,
776 ! CALIFORNIA INSTITUTE OF TECHNOLOGY, PASADENA, CALIFORNIA.
777 !.....RESTRICTIONS...
778 ! 1. THE MODEL EDDY DIFFUSIVITIES ARE BASED ON MONIN-OBUKHOV
779 ! SIMILARITY THEORY AND SO ARE ONLY APPLICABLE IN THE
780 ! SURFACE LAYER, A HEIGHT OF O(30M).
781 ! 2. ALL INPUT UNITS MUST BE CONSISTENT
782 ! 3. THE PHI FUNCTIONS USED TO CALCULATE THE FRICTION
783 ! VELOCITY U* AND THE POLLUTANT INTEGRALS ARE BASED
784 ! ON THE WORK OF BUSINGER ET AL.(1971).
785 ! 4. THE MOMENTUM AND POLLUTANT DIFFUSIVITIES ARE NOT
786 ! THE SAME FOR THE CASES L<0 AND L>0.
787 USE module_model_constants, only: karman
789 ! .. Scalar Arguments ..
790 REAL :: rmol, ustar, vgpart, z0, zr, aer_res
793 ! .. Array Arguments ..
796 ! .. Local Scalars ..
797 REAL :: ao, ar, polint, vk
800 ! .. Intrinsic Functions ..
803 ! Set the von Karman constant
806 ! Calculate the diffusion correction factor
807 ! SCPR23 is calculated as (Sc/Pr)**(2/3) using Sc= 1.15 and Pr= 1.0
810 ! DETERMINE THE STABILITY BASED ON THE CONDITIONS
815 if(abs(rmol) < 1.E-6 ) rmol = 0.
818 ar = ((1.0-9.0*zr*rmol)**(0.25)+0.001)**2
819 ao = ((1.0-9.0*z0*rmol)**(0.25)+0.001)**2
820 polint = 0.74*(alog((ar-1.0)/(ar+1.0))-alog((ao-1.0)/(ao+1.0)))
821 ELSE IF (rmol==0.) THEN
822 polint = 0.74*alog(zr/z0)
824 polint = 0.74*alog(zr/z0) + 4.7*rmol*(zr-z0)
827 ! CALCULATE THE Maximum DEPOSITION VELOCITY
830 depv(l) = ustar*vk/(2.0*scpr23(l)+polint)
832 vgpart = ustar*vk/polint
833 aer_res = polint/(karman*max(ustar,1.0e-4))
836 END SUBROUTINE depvel
838 SUBROUTINE dep_init(id,config_flags,numgas)
839 USE module_model_constants
841 USE module_state_description
842 TYPE (grid_config_rec_type) , INTENT (in) :: config_flags
844 ! .. Scalar Arguments ..
845 integer, intent(in) :: id, numgas
848 ! .. Local Scalars ..
850 INTEGER :: iland, iseason, l
854 REAL :: dat1(nlu,dep_seasons), dat2(nlu,dep_seasons), &
855 dat3(nlu,dep_seasons), dat4(nlu,dep_seasons), &
856 dat5(nlu,dep_seasons), dat6(nlu,dep_seasons), &
857 dat7(nlu,dep_seasons), dvj(numgas)
859 ! .. Make sure that the model is being run with a soil model. Otherwise,
860 ! iland will be zero in deppart, which will try to pull non-exisant
862 call nl_get_sf_surface_physics(id,l)
864 call wrf_error_fatal("ERROR: Cannot use dry deposition without using a soil model.")
867 ! .. Data Statements ..
868 ! RI for stomatal resistance
869 ! data ((ri(ILAND,ISEASON),ILAND=1,nlu),ISEASON=1,dep_seasons)/0.10E+11, &
870 DATA ((dat1(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+11, &
871 0.60E+02, 0.60E+02, 0.60E+02, 0.60E+02, 0.70E+02, 0.12E+03, &
872 0.12E+03, 0.12E+03, 0.12E+03, 0.70E+02, 0.13E+03, 0.70E+02, &
873 0.13E+03, 0.10E+03, 0.10E+11, 0.80E+02, 0.10E+03, 0.10E+11, &
874 0.80E+02, 0.10E+03, 0.10E+03, 0.10E+11, 0.10E+11, 0.10E+11, &
875 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
876 0.10E+11, 0.10E+11, 0.10E+11, 0.12E+03, 0.10E+11, 0.10E+11, &
877 0.70E+02, 0.25E+03, 0.50E+03, 0.10E+11, 0.10E+11, 0.50E+03, &
878 0.10E+11, 0.10E+11, 0.50E+03, 0.50E+03, 0.10E+11, 0.10E+11, &
879 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
880 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.12E+03, 0.10E+11, &
881 0.10E+11, 0.70E+02, 0.25E+03, 0.50E+03, 0.10E+11, 0.10E+11, &
882 0.50E+03, 0.10E+11, 0.10E+11, 0.50E+03, 0.50E+03, 0.10E+11, &
883 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
884 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
885 0.10E+11, 0.10E+11, 0.70E+02, 0.40E+03, 0.80E+03, 0.10E+11, &
886 0.10E+11, 0.80E+03, 0.10E+11, 0.10E+11, 0.80E+03, 0.80E+03, &
887 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.12E+03, 0.12E+03, &
888 0.12E+03, 0.12E+03, 0.14E+03, 0.24E+03, 0.24E+03, 0.24E+03, &
889 0.12E+03, 0.14E+03, 0.25E+03, 0.70E+02, 0.25E+03, 0.19E+03, &
890 0.10E+11, 0.16E+03, 0.19E+03, 0.10E+11, 0.16E+03, 0.19E+03, &
891 0.19E+03, 0.10E+11, 0.10E+11, 0.10E+11/
894 call wrf_debug(0, 'number of land use classifications not correct ')
895 CALL wrf_error_fatal ( "LAND USE CLASSIFICATIONS NOT 25")
897 IF (dep_seasons/=5) THEN
898 call wrf_debug(0, 'number of dep_seasons not correct ')
899 CALL wrf_error_fatal ( "DEP_SEASONS NOT 5")
902 ! SURFACE RESISTANCE DATA FOR DEPOSITION MODEL OF
903 ! M. L. WESELY, ATMOSPHERIC ENVIRONMENT 23 (1989) 1293-1304
905 ! Seasonal categories:
906 ! 1: midsummer with lush vegetation
907 ! 2: autumn with unharvested cropland
908 ! 3: late autumn with frost, no snow
909 ! 4: winter, snow on ground and subfreezing
910 ! 5: transitional spring with partially green short annuals
913 ! USGS type Wesely type
914 ! 1: Urban and built-up land 1
915 ! 2: Dryland cropland and pasture 2
916 ! 3: Irrigated cropland and pasture 2
917 ! 4: Mix. dry/irrg. cropland and pasture 2
918 ! 5: Cropland/grassland mosaic 2
919 ! 6: Cropland/woodland mosaic 4
922 ! 9: Mixed shrubland/grassland 3
923 ! 10: Savanna 3, always summer
924 ! 11: Deciduous broadleaf forest 4
925 ! 12: Deciduous needleleaf forest 5, autumn and winter modi
926 ! 13: Evergreen broadleaf forest 4, always summer
927 ! 14: Evergreen needleleaf forest 5
930 ! 17: Herbaceous wetland 9
931 ! 18: Wooded wetland 6
932 ! 19: Barren or sparsely vegetated 8
933 ! 20: Herbaceous Tundra 9
934 ! 21: Wooded Tundra 6
936 ! 23: Bare Ground Tundra 8
937 ! 24: Snow or Ice -, always winter
943 ! | seasonal category
947 ! RLU for outer surfaces in the upper canopy
948 DO iseason = 1, dep_seasons
950 ri(iland,iseason) = dat1(iland,iseason)
953 ! data ((rlu(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.10E+11, &
954 DATA ((dat2(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+11, &
955 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, &
956 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, &
957 0.20E+04, 0.20E+04, 0.10E+11, 0.25E+04, 0.20E+04, 0.10E+11, &
958 0.25E+04, 0.20E+04, 0.20E+04, 0.10E+11, 0.10E+11, 0.10E+11, &
959 0.10E+11, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, &
960 0.90E+04, 0.90E+04, 0.90E+04, 0.20E+04, 0.90E+04, 0.90E+04, &
961 0.20E+04, 0.40E+04, 0.80E+04, 0.10E+11, 0.90E+04, 0.80E+04, &
962 0.10E+11, 0.90E+04, 0.80E+04, 0.80E+04, 0.10E+11, 0.10E+11, &
963 0.10E+11, 0.10E+11, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, &
964 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, 0.20E+04, 0.90E+04, &
965 0.90E+04, 0.20E+04, 0.40E+04, 0.80E+04, 0.10E+11, 0.90E+04, &
966 0.80E+04, 0.10E+11, 0.90E+04, 0.80E+04, 0.80E+04, 0.10E+11, &
967 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
968 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
969 0.10E+11, 0.10E+11, 0.20E+04, 0.60E+04, 0.90E+04, 0.10E+11, &
970 0.90E+04, 0.90E+04, 0.10E+11, 0.90E+04, 0.90E+04, 0.90E+04, &
971 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.40E+04, 0.40E+04, &
972 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, &
973 0.20E+04, 0.40E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.30E+04, &
974 0.10E+11, 0.40E+04, 0.30E+04, 0.10E+11, 0.40E+04, 0.30E+04, &
975 0.30E+04, 0.10E+11, 0.10E+11, 0.10E+11/
976 DO iseason = 1, dep_seasons
978 rlu(iland,iseason) = dat2(iland,iseason)
981 ! RAC for transfer that depends on canopy height and density
982 ! data ((rac(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.10E+03, &
983 DATA ((dat3(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+03, &
984 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+04, 0.10E+03, &
985 0.10E+03, 0.10E+03, 0.10E+03, 0.20E+04, 0.20E+04, 0.20E+04, &
986 0.20E+04, 0.20E+04, 0.00E+00, 0.30E+03, 0.20E+04, 0.00E+00, &
987 0.30E+03, 0.20E+04, 0.20E+04, 0.00E+00, 0.00E+00, 0.00E+00, &
988 0.10E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+04, &
989 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, 0.15E+04, 0.20E+04, &
990 0.20E+04, 0.20E+04, 0.17E+04, 0.00E+00, 0.20E+03, 0.17E+04, &
991 0.00E+00, 0.20E+03, 0.17E+04, 0.17E+04, 0.00E+00, 0.00E+00, &
992 0.00E+00, 0.10E+03, 0.10E+02, 0.10E+02, 0.10E+02, 0.10E+02, &
993 0.10E+04, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+04, &
994 0.20E+04, 0.20E+04, 0.20E+04, 0.15E+04, 0.00E+00, 0.10E+03, &
995 0.15E+04, 0.00E+00, 0.10E+03, 0.15E+04, 0.15E+04, 0.00E+00, &
996 0.00E+00, 0.00E+00, 0.10E+03, 0.10E+02, 0.10E+02, 0.10E+02, &
997 0.10E+02, 0.10E+04, 0.10E+02, 0.10E+02, 0.10E+02, 0.10E+02, &
998 0.10E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.15E+04, 0.00E+00, &
999 0.50E+02, 0.15E+04, 0.00E+00, 0.50E+02, 0.15E+04, 0.15E+04, &
1000 0.00E+00, 0.00E+00, 0.00E+00, 0.10E+03, 0.50E+02, 0.50E+02, &
1001 0.50E+02, 0.50E+02, 0.12E+04, 0.80E+02, 0.80E+02, 0.80E+02, &
1002 0.10E+03, 0.12E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.15E+04, &
1003 0.00E+00, 0.20E+03, 0.15E+04, 0.00E+00, 0.20E+03, 0.15E+04, &
1004 0.15E+04, 0.00E+00, 0.00E+00, 0.00E+00/
1005 DO iseason = 1, dep_seasons
1007 rac(iland,iseason) = dat3(iland,iseason)
1010 ! RGSS for ground surface SO2
1011 ! data ((rgss(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.40E+03, &
1012 DATA ((dat4(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.40E+03, &
1013 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.50E+03, 0.35E+03, &
1014 0.35E+03, 0.35E+03, 0.35E+03, 0.50E+03, 0.50E+03, 0.50E+03, &
1015 0.50E+03, 0.10E+03, 0.10E+01, 0.10E+01, 0.10E+03, 0.10E+04, &
1016 0.10E+01, 0.10E+03, 0.10E+03, 0.10E+04, 0.10E+03, 0.10E+04, &
1017 0.40E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.50E+03, &
1018 0.35E+03, 0.35E+03, 0.35E+03, 0.35E+03, 0.50E+03, 0.50E+03, &
1019 0.50E+03, 0.50E+03, 0.10E+03, 0.10E+01, 0.10E+01, 0.10E+03, &
1020 0.10E+04, 0.10E+01, 0.10E+03, 0.10E+03, 0.10E+04, 0.10E+03, &
1021 0.10E+04, 0.40E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, &
1022 0.50E+03, 0.35E+03, 0.35E+03, 0.35E+03, 0.35E+03, 0.50E+03, &
1023 0.50E+03, 0.50E+03, 0.50E+03, 0.20E+03, 0.10E+01, 0.10E+01, &
1024 0.20E+03, 0.10E+04, 0.10E+01, 0.20E+03, 0.20E+03, 0.10E+04, &
1025 0.10E+03, 0.10E+04, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, &
1026 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, &
1027 0.10E+03, 0.10E+03, 0.50E+03, 0.10E+03, 0.10E+03, 0.10E+01, &
1028 0.10E+03, 0.10E+03, 0.10E+04, 0.10E+03, 0.10E+03, 0.10E+03, &
1029 0.10E+04, 0.10E+03, 0.10E+04, 0.50E+03, 0.15E+03, 0.15E+03, &
1030 0.15E+03, 0.15E+03, 0.50E+03, 0.35E+03, 0.35E+03, 0.35E+03, &
1031 0.35E+03, 0.50E+03, 0.50E+03, 0.50E+03, 0.50E+03, 0.20E+03, &
1032 0.10E+01, 0.10E+01, 0.20E+03, 0.10E+04, 0.10E+01, 0.20E+03, &
1033 0.20E+03, 0.10E+04, 0.10E+03, 0.10E+04/
1034 DO iseason = 1, dep_seasons
1036 rgss(iland,iseason) = dat4(iland,iseason)
1039 ! RGSO for ground surface O3
1040 ! data ((rgso(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.30E+03, &
1041 DATA ((dat5(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.30E+03, &
1042 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.20E+03, 0.20E+03, &
1043 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, &
1044 0.20E+03, 0.30E+03, 0.20E+04, 0.10E+04, 0.30E+03, 0.40E+03, &
1045 0.10E+04, 0.30E+03, 0.30E+03, 0.40E+03, 0.35E+04, 0.40E+03, &
1046 0.30E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.20E+03, &
1047 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, &
1048 0.20E+03, 0.20E+03, 0.30E+03, 0.20E+04, 0.80E+03, 0.30E+03, &
1049 0.40E+03, 0.80E+03, 0.30E+03, 0.30E+03, 0.40E+03, 0.35E+04, &
1050 0.40E+03, 0.30E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, &
1051 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, &
1052 0.20E+03, 0.20E+03, 0.20E+03, 0.30E+03, 0.20E+04, 0.10E+04, &
1053 0.30E+03, 0.40E+03, 0.10E+04, 0.30E+03, 0.30E+03, 0.40E+03, &
1054 0.35E+04, 0.40E+03, 0.60E+03, 0.35E+04, 0.35E+04, 0.35E+04, &
1055 0.35E+04, 0.35E+04, 0.35E+04, 0.35E+04, 0.35E+04, 0.35E+04, &
1056 0.35E+04, 0.35E+04, 0.20E+03, 0.35E+04, 0.35E+04, 0.20E+04, &
1057 0.35E+04, 0.35E+04, 0.40E+03, 0.35E+04, 0.35E+04, 0.35E+04, &
1058 0.40E+03, 0.35E+04, 0.40E+03, 0.30E+03, 0.15E+03, 0.15E+03, &
1059 0.15E+03, 0.15E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, &
1060 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.30E+03, &
1061 0.20E+04, 0.10E+04, 0.30E+03, 0.40E+03, 0.10E+04, 0.30E+03, &
1062 0.30E+03, 0.40E+03, 0.35E+04, 0.40E+03/
1063 DO iseason = 1, dep_seasons
1065 rgso(iland,iseason) = dat5(iland,iseason)
1068 ! RCLS for exposed surfaces in the lower canopy SO2
1069 ! data ((rcls(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.10E+11, &
1070 DATA ((dat6(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+11, &
1071 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, &
1072 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, &
1073 0.20E+04, 0.20E+04, 0.10E+11, 0.25E+04, 0.20E+04, 0.10E+11, &
1074 0.25E+04, 0.20E+04, 0.20E+04, 0.10E+11, 0.10E+11, 0.10E+11, &
1075 0.10E+11, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, &
1076 0.90E+04, 0.90E+04, 0.90E+04, 0.20E+04, 0.90E+04, 0.90E+04, &
1077 0.20E+04, 0.20E+04, 0.40E+04, 0.10E+11, 0.90E+04, 0.40E+04, &
1078 0.10E+11, 0.90E+04, 0.40E+04, 0.40E+04, 0.10E+11, 0.10E+11, &
1079 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
1080 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, 0.20E+04, 0.90E+04, &
1081 0.90E+04, 0.20E+04, 0.30E+04, 0.60E+04, 0.10E+11, 0.90E+04, &
1082 0.60E+04, 0.10E+11, 0.90E+04, 0.60E+04, 0.60E+04, 0.10E+11, &
1083 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
1084 0.10E+11, 0.90E+04, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
1085 0.90E+04, 0.90E+04, 0.20E+04, 0.20E+03, 0.40E+03, 0.10E+11, &
1086 0.90E+04, 0.40E+03, 0.10E+11, 0.90E+04, 0.40E+03, 0.40E+03, &
1087 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.40E+04, 0.40E+04, &
1088 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, &
1089 0.20E+04, 0.40E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.30E+04, &
1090 0.10E+11, 0.40E+04, 0.30E+04, 0.10E+11, 0.40E+04, 0.30E+04, &
1091 0.30E+04, 0.10E+11, 0.10E+11, 0.10E+11/
1092 DO iseason = 1, dep_seasons
1094 rcls(iland,iseason) = dat6(iland,iseason)
1097 ! RCLO for exposed surfaces in the lower canopy O3
1098 ! data ((rclo(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.10E+11, &
1099 DATA ((dat7(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+11, &
1100 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, &
1101 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, &
1102 0.10E+04, 0.10E+04, 0.10E+11, 0.10E+04, 0.10E+04, 0.10E+11, &
1103 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+11, 0.10E+11, 0.10E+11, &
1104 0.10E+11, 0.40E+03, 0.40E+03, 0.40E+03, 0.40E+03, 0.40E+03, &
1105 0.40E+03, 0.40E+03, 0.40E+03, 0.10E+04, 0.40E+03, 0.40E+03, &
1106 0.10E+04, 0.10E+04, 0.60E+03, 0.10E+11, 0.40E+03, 0.60E+03, &
1107 0.10E+11, 0.40E+03, 0.60E+03, 0.60E+03, 0.10E+11, 0.10E+11, &
1108 0.10E+11, 0.10E+11, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, &
1109 0.40E+03, 0.40E+03, 0.40E+03, 0.40E+03, 0.10E+04, 0.40E+03, &
1110 0.40E+03, 0.10E+04, 0.10E+04, 0.60E+03, 0.10E+11, 0.80E+03, &
1111 0.60E+03, 0.10E+11, 0.80E+03, 0.60E+03, 0.60E+03, 0.10E+11, &
1112 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+04, 0.10E+04, 0.10E+04, &
1113 0.10E+04, 0.40E+03, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, &
1114 0.40E+03, 0.40E+03, 0.10E+04, 0.15E+04, 0.60E+03, 0.10E+11, &
1115 0.80E+03, 0.60E+03, 0.10E+11, 0.80E+03, 0.60E+03, 0.60E+03, &
1116 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+04, 0.10E+04, &
1117 0.10E+04, 0.10E+04, 0.50E+03, 0.50E+03, 0.50E+03, 0.50E+03, &
1118 0.10E+04, 0.50E+03, 0.15E+04, 0.10E+04, 0.15E+04, 0.70E+03, &
1119 0.10E+11, 0.60E+03, 0.70E+03, 0.10E+11, 0.60E+03, 0.70E+03, &
1120 0.70E+03, 0.10E+11, 0.10E+11, 0.10E+11/
1121 DO iseason = 1, dep_seasons
1123 rclo(iland,iseason) = dat7(iland,iseason)
1134 ! HENRY''S LAW COEFFICIENTS
1135 ! Effective Henry''s law coefficient at pH 7
1136 ! [KH298]=mole/(l atm)
1138 if (config_flags%chem_opt /= CBM4_KPP) then
1140 hstar(p_no2) = 6.40E-3
1141 hstar(p_no) = 1.90E-3
1142 hstar(p_pan) = 2.97E+0
1143 hstar(p_o3) = 1.13E-2
1144 hstar(p_hcho) = 2.97E+3
1145 hstar(p_aco3) = 1.14E+1
1146 hstar(p_tpan) = 2.97E+0
1147 hstar(p_hono) = 3.47E+5
1148 hstar(p_no3) = 1.50E+1
1149 hstar(p_hno4) = 2.00E+13
1150 hstar(p_h2o2) = 7.45E+4
1151 hstar(p_co) = 8.20E-3
1152 hstar(p_ald) = 1.14E+1
1153 hstar(p_op1) = 2.21E+2
1154 hstar(p_op2) = 1.68E+6
1155 hstar(p_paa) = 4.73E+2
1156 hstar(p_ket) = 3.30E+1
1157 hstar(p_gly) = 1.40E+6
1158 hstar(p_mgly) = 3.71E+3
1159 hstar(p_dcb) = 1.40E+6
1160 hstar(p_onit) = 1.13E+0
1161 hstar(p_so2) = 2.53E+5
1162 hstar(p_eth) = 2.00E-3
1163 hstar(p_hc3) = 1.42E-3
1164 hstar(p_hc5) = 1.13E-3
1165 hstar(p_hc8) = 1.42E-3
1166 hstar(p_olt) = 4.76E-3
1167 hstar(p_oli) = 1.35E-3
1168 hstar(p_tol) = 1.51E-1
1169 hstar(p_csl) = 4.00E+5
1170 hstar(p_xyl) = 1.45E-1
1171 hstar(p_iso) = 4.76E-3
1172 hstar(p_hno3) = 2.69E+13
1173 hstar(p_ora1) = 9.85E+6
1174 hstar(p_ora2) = 9.63E+5
1175 hstar(p_nh3) = 1.04E+4
1176 hstar(p_n2o5) = 1.00E+10
1177 if(p_ol2.gt.1) hstar(p_ol2) = 4.67E-3
1178 if(p_par.gt.1) hstar(p_par) = 1.13E-3 !wig, 1-May-2007: for CBMZ
1180 hstar(p_ch4) = 1.50E-3
1186 hstar(p_co2) = 1.86E-1
1192 ! FOLLOWING FOR RACM
1195 HSTAR(p_ETE )=4.67E-3
1196 HSTAR(p_API )=4.76E-3
1197 HSTAR(p_LIM )=4.76E-3
1198 HSTAR(p_DIEN)=4.76E-3
1199 HSTAR(p_MACR)=1.14E+1
1200 HSTAR(p_UDD )=1.40E+6
1201 HSTAR(p_HKET)=7.80E+3
1224 ! -DH/R (for temperature correction)
1264 if(p_ol2.gt.1)dhr(p_ol2) = 0.
1265 if(p_par.gt.1)dhr(p_par) = 0. !wig, 1-May-2007: for CBMZ
1266 ! REACTIVITY FACTORS
1306 if(p_ol2.gt.1)f0(p_ol2) = 0.
1307 if(p_par.gt.1)f0(p_par) = 0. !wig, 1-May-2007: for CBMZ
1308 ! DIFFUSION COEFFICIENTS
1309 ! [DV]=cm2/s (assumed: 1/SQRT(molar mass) when not known)
1350 if(p_ol2.gt.1)dvj(p_ol2) = 0.189
1351 if(p_par.gt.1)dvj(p_par) = 0.118 !wig, 1-May-2007: for CBMZ
1353 hstar4(l) = hstar(l) ! preliminary
1354 ! Correction of diff. coeff
1355 dvj(l) = dvj(l)*(293.15/298.15)**1.75
1356 sc = 0.15/dvj(l) ! Schmidt Number at 20°C
1357 dratio(l) = 0.242/dvj(l) ! ! of water vapor and gas at
1358 ! Ratio of diffusion coeffi
1359 scpr23(l) = (sc/0.72)**(2./3.) ! (Schmidt # / Prandtl #)**
1363 hstar(p_no2) = 6.40E-3
1364 hstar(p_no) = 1.90E-3
1365 hstar(p_pan) = 2.97E+0
1366 hstar(p_o3) = 1.13E-2
1367 hstar(p_hcho) = 2.97E+3
1368 hstar(p_hono) = 3.47E+5
1369 hstar(p_no3) = 1.50E+1
1370 hstar(p_h2o2) = 7.45E+4
1371 hstar(p_co) = 8.20E-3
1372 hstar(p_ald2) = 1.14E+1
1373 hstar(p_onit) = 1.13E+0
1374 hstar(p_so2) = 2.53E+5
1375 hstar(p_eth) = 2.00E-3
1376 hstar(p_ole) = 3.05E-3 !(average of oli and olt)
1377 hstar(p_tol) = 1.51E-1
1378 hstar(p_cres) = 4.00E+5
1379 hstar(p_xyl) = 1.45E-1
1380 hstar(p_iso) = 4.76E-3
1381 hstar(p_hno3) = 2.69E+13
1382 hstar(p_nh3) = 1.04E+4
1383 hstar(p_n2o5) = 1.00E+10
1384 hstar(p_par) = 1.13E-3
1386 ! -DH/R (for temperature correction)
1411 ! REACTIVITY FACTORS
1436 ! DIFFUSION COEFFICIENTS
1437 ! [DV]=cm2/s (assumed: 1/SQRT(molar mass) when not known)
1462 hstar4(l) = hstar(l) ! preliminary
1463 ! Correction of diff. coeff
1464 dvj(l) = dvj(l)*(293.15/298.15)**1.75
1465 sc = 0.15/dvj(l) ! Schmidt Number at 20°C
1466 dratio(l) = 0.242/dvj(l) ! ! of water vapor and gas at
1467 ! Ratio of diffusion coeffi
1468 scpr23(l) = (sc/0.72)**(2./3.) ! (Schmidt # / Prandtl #)**
1475 ! DATA FOR AEROSOL PARTICLE DEPOSITION FOR THE MODEL OF
1476 ! J. W. ERISMAN, A. VAN PUL AND P. WYERS
1477 ! ATMOSPHERIC ENVIRONMENT 28 (1994), 2595-2607
1479 ! vd = (u* / k) * CORRECTION FACTORS
1481 ! CONSTANT K FOR LANDUSE TYPES:
1482 ! urban and built-up land
1484 ! dryland cropland and pasture
1486 ! irrigated cropland and pasture
1488 ! mixed dryland/irrigated cropland and past
1490 ! cropland/grassland mosaic
1492 ! cropland/woodland mosaic
1498 ! mixed shrubland/grassland
1502 ! deciduous broadleaf forest
1504 ! deciduous needleleaf forest
1506 ! evergreen broadleaf forest
1508 ! evergreen needleleaf forest
1514 ! herbaceous wetland
1518 ! barren or sparsely vegetated
1526 ! bare ground tundra
1532 ! Erisman et al. (1994) give
1533 ! k = 500 for low vegetation and k = 100 for forests.
1535 ! For desert k = 500 is taken according to measurements
1537 ! J. Fontan, A. Lopez, E. Lamaud and A. Druilhet (1997)
1538 ! Vertical Flux Measurements of the Submicronic Aerosol Particles
1539 ! and Parametrisation of the Dry Deposition Velocity
1540 ! in: Biosphere-Atmosphere Exchange of Pollutants
1541 ! and Trace Substances
1542 ! Editor: S. Slanina. Springer-Verlag Berlin, Heidelberg, 1997
1545 ! For coniferous forest the Erisman value of k = 100 is taken.
1546 ! Measurements of Erisman et al. (1997) in a coniferous forest
1547 ! in the Netherlands, lead to values of k between 20 and 38
1548 ! (Atmospheric Environment 31 (1997), 321-332).
1549 ! However, these high values of vd may be reached during
1550 ! instable cases. The eddy correlation measurements
1551 ! of Gallagher et al. (1997) made during the same experiment
1552 ! show for stable cases (L>0) values of k between 200 and 250
1553 ! at minimum (Atmospheric Environment 31 (1997), 359-373).
1554 ! Fontan et al. (1997) found k = 250 in a forest
1555 ! of maritime pine in southwestern France.
1557 ! For gras, model calculations of Davidson et al. support
1559 ! C. I. Davidson, J. M. Miller and M. A. Pleskov
1560 ! The Influence of Surface Structure on Predicted Particles
1561 ! Dry Deposition to Natural Gras Canopies
1562 ! Water, Air, and Soil Pollution 18 (1982) 25-43
1564 ! Snow covered surface: The experiment of Ibrahim et al. (1983)
1565 ! gives k = 436 for 0.7 um diameter particles.
1566 ! The deposition velocity of Milford and Davidson (1987)
1567 ! gives k = 154 for continental sulfate aerosol.
1568 ! M. Ibrahim, L. A. Barrie and F. Fanaki
1569 ! Atmospheric Environment 17 (1983), 781-788
1571 ! J. B. Milford and C. I. Davidson
1572 ! The Sizes of Particulate Sulfate and Nitrate in the Atmosphere
1574 ! JAPCA 37 (1987), 125-134
1576 ! WRITE (0,*) ' return from rcread '
1577 ! *********************************************************
1579 ! Simplified landuse scheme for deposition and biogenic emission
1581 ! (ISWATER and ISICE are already defined elsewhere,
1582 ! therefore water and ice are not considered here)
1584 ! 1 urban or bare soil
1587 ! 4 deciduous forest
1588 ! 5 coniferous and mixed forest
1589 ! 6 other natural landuse categories
1592 IF (mminlu=='OLD ') THEN
1607 IF (mminlu=='USGS') THEN
1634 IF (mminlu=='SiB ') THEN
1654 END SUBROUTINE dep_init
1655 END MODULE module_dep_simple