module_dep_simple.F

References to this file elsewhere.
1     MODULE module_dep_simple
2       IMPLICIT NONE
3 !
4 ! many of these parameters will depend on the RADM mechanism!
5 ! if you change it, lets talk about it and get it done!!!
6 !
7 
8       INTEGER, PARAMETER :: dep_seasons = 5 ,   nlu = 25,                   &
9         nseason = 1, nseasons = 2
10 !
11 ! following currently hardwired to USGS
12 !
13       INTEGER, PARAMETER :: isice_temp=24,iswater_temp=16
14       character, parameter :: mminlu='USGS'
15 !
16       INTEGER :: ixxxlu(nlu)
17       REAL :: kpart(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)
21 !
22 ! NO MORE THAN 1000 SPECIES FOR DEPOSITION
23 !
24       REAL, DIMENSION (1:1000) :: dratio,hstar,hstar4,f0,dhr,scpr23
25 ! .. Default Accessibility ..
26     PUBLIC
27             logical, allocatable :: is_aerosol(:) ! true if field is aerosol (any phase)
28 ! ..
29     CONTAINS
30        subroutine wesely_driver(id,ktau,dtstep,                           &
31                config_flags,                                              &
32                gmt,julday,t_phy,moist,p8w,t8w,                            &
33                p_phy,chem,rho_phy,dz8w,ddvel,aer_res,                     &
34                ivgtyp,tsk,gsw,vegfra,pbl,rmol,ust,znt,xlat,xlong,z,z_at_w,&
35                numgas,                                                    &
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 
41   USE module_configure
42   USE module_state_description                       
43   USE module_data_sorgam
44 ! USE module_data_radm2                               
45 ! USE module_data_racm                               
46           
47    INTEGER,      INTENT(IN   ) :: id,julday,                              &
48                                   numgas,                                 &
49                                   ids,ide, jds,jde, kds,kde,              &
50                                   ims,ime, jms,jme, kms,kme,              &
51                                   its,ite, jts,jte, kts,kte     
52    INTEGER,      INTENT(IN   ) ::                                         &
53                                   ktau            
54       REAL,      INTENT(IN   ) ::                                         &
55                              dtstep,gmt
56 !
57 ! advected moisture variables
58 !
59    REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_moist ),               &
60          INTENT(IN ) ::                                   moist  
61 !
62 ! advected chemical species
63 !
64    REAL, DIMENSION( ims:ime, kms:kme, jms:jme, num_chem ),                &
65          INTENT(INOUT ) ::                                   chem
66 !
67 ! deposition velocities
68 !
69    REAL, DIMENSION( its:ite, jts:jte, num_chem ),                         &
70          INTENT(INOUT ) ::   ddvel                     
71 !
72 ! input from met model
73 !
74 
75    REAL,  DIMENSION( ims:ime , kms:kme , jms:jme )         ,              &
76           INTENT(IN   ) ::                                                &
77                                                       t_phy,              &
78                                                       p_phy,              &
79                                                       dz8w,               &
80                                                       z    ,              &
81                                               t8w,p8w,z_at_w ,            &
82                                                     rho_phy
83    INTEGER,DIMENSION( ims:ime , jms:jme )                  ,              &
84           INTENT(IN   ) ::                                                &
85                                                      ivgtyp
86    REAL,  DIMENSION( ims:ime , jms:jme )                   ,              &
87           INTENT(INOUT   ) ::                                                &
88                                                      tsk,                 &
89                                                      gsw,                 &
90                                                   vegfra,                 &
91                                                      pbl,                 &
92                                                      rmol,                &
93                                                      ust,                 &
94                                                      xlat,                &
95                                                      xlong,               &
96                                                      znt
97 
98 !--- deposition and emissions stuff
99 ! ..
100 ! .. Local Scalars ..
101       REAL ::  clwchem,dvfog,dvpart,pa,rad,rhchem,ta,ustar,vegfrac,z1,zntt
102       INTEGER :: i, iland, iprt, iseason, j, jce, jcs, n, nr, ipr,jpr,nvr
103       LOGICAL :: highnh3, rainflag, vegflag, wetflag
104 ! ..
105 ! .. Local Arrays ..  
106       REAL :: p(kts:kte), srfres(numgas),ddvel0d(numgas)
107 !
108 ! necessary for aerosols (module dependent)         
109 !  
110       real :: aer_res(its:ite,jts:jte),rcx(numgas)
111                                                      
112    TYPE(grid_config_rec_type),  INTENT(IN   )    :: config_flags
113           
114                                                      
115 ! ..                                                 
116 ! .. Intrinsic Functions ..                       
117       INTRINSIC max, min                             
118 ! ..                                                 
119       CALL wrf_debug(15,'in dry_dep_wesely')
120       iseason = 1
121       if(julday.lt.90.or.julday.gt.270)then
122         iseason=2
123         CALL wrf_debug(15,'setting iseason to 2')
124       endif
125       do 100 j=jts,jte 
126       do 100 i=its,ite 
127       iprt=0
128       iland = ivgtyp(i,j)
129       ta = tsk(i,j)
130       rad = gsw(i,j)
131       vegfrac = vegfra(i,j)
132       pa = .01*p_phy(i,kts,j)
133       clwchem = moist(i,kts,j,p_qc)
134       ustar = ust(i,j)
135       zntt = znt(i,j)
136 
137       z1 = z_at_w(i,kts+1,j)-z_at_w(i,kts,j)
138 
139 !     Set logical default values
140       rainflag = .FALSE.
141       wetflag = .FALSE.
142       highnh3 = .FALSE.
143       if(p_qr.gt.1)then
144          if(moist(i,kts,j,p_qr).gt.0.)rainflag = .true.
145       endif
146       rhchem = MIN( 100.,100. * moist(i,kts,j,p_qv) / &
147                (3.80*exp(17.27*(t_phy(i,kts,j)-273.)/(t_phy(i,kts,j)-36.))/pa))
148       rhchem = MAX(5.,RHCHEM)
149       if (rhchem >= 95.) wetflag = .true.
150 
151       if(chem(i,kts,j,p_nh3).gt.2.*chem(i,kts,j,p_so2))highnh3 = .true.
152 
153 !--- deposition
154       
155 !     if(snowc(i,j).gt.0.)iseason=4
156       CALL rc(rcx,ta,rad,rhchem,iland,iseason,numgas,      &
157         wetflag,rainflag,highnh3,iprt,p_o3,p_so2,p_nh3)
158       srfres=0.
159       DO n = 1, numgas-2
160         srfres(n) = rcx(n)
161       END DO
162       CALL deppart(rmol(i,j),ustar,rhchem,clwchem,iland,dvpart,dvfog)
163       ddvel0d=0.
164       aer_res(i,j)=0.
165       CALL landusevg(ddvel0d,ustar,rmol(i,j),zntt,z1,dvpart,iland,        &
166            numgas,srfres,aer_res(i,j),p_sulf)
167 
168 !wig: CMBZ does not have HO and HO2 last so need to copy all species
169 !      ddvel(i,j,1:numgas-2)=ddvel0d(1:numgas-2)
170       ddvel(i,j,1:numgas)=ddvel0d(1:numgas)
171 
172 100   continue
173 !
174 ! For the additional CBMZ species, assign similar RADM counter parts for
175 ! now. Short lived species get a zero velocity since dry dep should be
176 ! unimportant.  **ALSO**, treat p_sulf as h2so4 vapor, not aerosol sulfate
177 !
178       if ( (config_flags%chem_opt == CBMZ          ) .or.   &
179            (config_flags%chem_opt == CBMZ_BB       ) .or.   &
180            (config_flags%chem_opt == CBMZ_MOSAIC_4BIN_AQ) .or.   &
181            (config_flags%chem_opt == CBMZ_MOSAIC_8BIN_AQ) .or.   &
182            (config_flags%chem_opt == CBMZ_MOSAIC_4BIN) .or.   &
183            (config_flags%chem_opt == CBMZ_MOSAIC_8BIN) ) then
184          do j=jts,jte
185             do i=its,ite
186                ddvel(i,j,p_sulf)        = ddvel(i,j,p_hno3)
187                ddvel(i,j,p_hcl)         = ddvel(i,j,p_hno3)
188                ddvel(i,j,p_hcl)         = ddvel(i,j,p_hno3)
189                ddvel(i,j,p_ch3o2)       = 0
190                ddvel(i,j,p_ethp)        = 0
191                ddvel(i,j,p_ch3oh)       = ddvel(i,j,p_hcho)
192                ddvel(i,j,p_c2h5oh)      = ddvel(i,j,p_hcho)
193 !wig, 1-May-2007 (added par to tables)    ddvel(i,j,p_par)         = ddvel(i,j,p_hc5)
194                ddvel(i,j,p_to2)         = 0
195                ddvel(i,j,p_cro)         = 0
196                ddvel(i,j,p_open)        = ddvel(i,j,p_xyl)
197                ddvel(i,j,p_op3)         = ddvel(i,j,p_op2)
198                ddvel(i,j,p_c2o3)        = 0
199                ddvel(i,j,p_ro2)         = 0
200                ddvel(i,j,p_ano2)        = 0
201                ddvel(i,j,p_nap)         = 0
202                ddvel(i,j,p_xo2)         = 0
203                ddvel(i,j,p_xpar)        = 0
204                ddvel(i,j,p_isoprd)      = 0
205                ddvel(i,j,p_isopp)       = 0
206                ddvel(i,j,p_isopn)       = 0
207                ddvel(i,j,p_isopo2)      = 0
208                if( config_flags%chem_opt == CBMZ ) then
209                   ddvel(i,j,p_dms)         = 0
210                   ddvel(i,j,p_msa)         = ddvel(i,j,p_hno3)
211                   ddvel(i,j,p_dmso)        = 0
212                   ddvel(i,j,p_dmso2)       = 0
213                   ddvel(i,j,p_ch3so2h)     = 0
214                   ddvel(i,j,p_ch3sch2oo)   = 0
215                   ddvel(i,j,p_ch3so2)      = 0
216                   ddvel(i,j,p_ch3so3)      = 0
217                   ddvel(i,j,p_ch3so2oo)    = 0
218                   ddvel(i,j,p_ch3so2ch2oo) = 0
219                   ddvel(i,j,p_mtf)         = 0
220                end if
221             end do
222          end do
223       end if
224 
225 END SUBROUTINE wesely_driver
226 
227 ! **********************************************************************
228 ! **************************  SUBROUTINE RC  ***************************
229 ! **********************************************************************
230       SUBROUTINE rc(rcx,t,rad,rh,iland,iseason,numgas,             &
231           wetflag,rainflag,highnh3,iprt,p_o3,p_so2,p_nh3)
232 !     THIS SUBROUTINE CALCULATES SURFACE RESISTENCES ACCORDING
233 !     TO THE MODEL OF
234 !     M. L. WESELY,
235 !     ATMOSPHERIC ENVIRONMENT 23 (1989), 1293-1304
236 !     WITH SOME ADDITIONS ACCORDING TO
237 !     J. W. ERISMAN, A. VAN PUL, AND P. WYERS,
238 !     ATMOSPHERIC ENVIRONMENT 28 (1994), 2595-2607
239 !     WRITTEN BY  WINFRIED SEIDL, APRIL 1997
240 !     MODYFIED BY WINFRIED SEIDL, MARCH 2000
241 !                    FOR MM5 VERSION 3
242 !----------------------------------------------------------------------
243 ! .. Scalar Arguments ..
244         REAL :: rad, rh, t
245         INTEGER :: iland, iseason, numgas
246         LOGICAL :: highnh3, rainflag, wetflag
247         real :: rcx(numgas)
248 ! ..
249         INTEGER :: iprt,p_o3,p_so2,p_nh3
250 ! ..
251 ! .. Local Scalars ..
252         REAL :: rclx, rdc, resice, rgsx, rluo1, rluo2, rlux, rmx, rs, rsmx, &
253           tc, rdtheta, z
254         INTEGER :: n
255 ! ..
256 ! .. Local Arrays ..
257         REAL :: hstary(numgas)
258 ! ..
259 ! .. Intrinsic Functions ..
260         INTRINSIC exp
261 ! ..
262         DO n = 1, numgas
263           rcx(n) = 1.
264         END DO
265 
266         tc = t - 273.15
267         rdtheta = 0.
268 
269         z = 200./(rad+0.1)
270 
271 !!!  HARDWIRE VALUES FOR TESTING
272 !       z=0.4727409
273 !       tc=22.76083
274 !       t=tc+273.15
275 !       rad = 412.8426
276 !       rainflag=.false.
277 !       wetflag=.false.
278 
279         IF ((tc<=0.) .OR. (tc>=40.)) THEN
280           rs = 9999.
281         ELSE
282           rs = ri(iland,iseason)*(1+z*z)*(400./(tc*(40.-tc)))
283         END IF
284         rdc = 100*(1.+1000./(rad+10))/(1+1000.*rdtheta)
285         rluo1 = 1./(1./3000.+1./3./rlu(iland,iseason))
286         rluo2 = 1./(1./1000.+1./3./rlu(iland,iseason))
287         resice = 1000.*exp(-tc-4.)
288 
289         DO n = 1, numgas
290           IF (hstar(n)==0.) GO TO 10
291           hstary(n) = hstar(n)*exp(dhr(n)*(1./t-1./298.))
292           rmx = 1./(hstary(n)/3000.+100.*f0(n))
293           rsmx = rs*dratio(n) + rmx
294           rclx = 1./(hstary(n)/1.E+5/rcls(iland,iseason)+f0(n)/rclo(iland, &
295             iseason)) + resice
296           rgsx = 1./(hstary(n)/1.E+5/rgss(iland,iseason)+f0(n)/rgso(iland, &
297             iseason)) + resice
298           rlux = rlu(iland,iseason)/(1.E-5*hstary(n)+f0(n)) + resice
299           IF (wetflag) THEN
300             rlux = 1./(1./3./rlu(iland,iseason)+1.E-7*hstary(n)+f0(n)/rluo1)
301           END IF
302           IF (rainflag) THEN
303             rlux = 1./(1./3./rlu(iland,iseason)+1.E-7*hstary(n)+f0(n)/rluo2)
304           END IF
305           rcx(n) = 1./(1./rsmx+1./rlux+1./(rdc+rclx)+1./(rac(iland, &
306             iseason)+rgsx))
307           IF (rcx(n)<1.) rcx(n) = 1.
308 10      END DO
309 
310 !     SPECIAL TREATMENT FOR OZONE
311         hstary(p_o3) = hstar(p_o3)*exp(dhr(p_o3)*(1./t-1./298.))
312         rmx = 1./(hstary(p_o3)/3000.+100.*f0(p_o3))
313         rsmx = rs*dratio(p_o3) + rmx
314         rlux = rlu(iland,iseason)/(1.E-5*hstary(p_o3)+f0(p_o3)) + resice
315         rclx = rclo(iland,iseason) + resice
316         rgsx = rgso(iland,iseason) + resice
317         IF (wetflag) rlux = rluo1
318         IF (rainflag) rlux = rluo2
319         rcx(p_o3) = 1./(1./rsmx+1./rlux+1./(rdc+rclx)+1./(rac(iland, &
320           iseason)+rgsx))
321         IF (rcx(p_o3)<1.) rcx(p_o3) = 1.
322 
323 !     SPECIAL TREATMENT FOR SO2 (Wesely)
324 !       HSTARY(P_SO2)=HSTAR(P_SO2)*EXP(DHR(P_SO2)*(1./T-1./298.))
325 !       RMX=1./(HSTARY(P_SO2)/3000.+100.*F0(P_SO2))
326 !       RSMX=RS*DRATIO(P_SO2)+RMX
327 !       RLUX=RLU(ILAND,ISEASON)/(1.E-5*HSTARY(P_SO2)+F0(P_SO2))
328 !    &       +RESICE
329 !       RCLX=RCLS(ILAND,ISEASON)+RESICE
330 !       RGSX=RGSS(ILAND,ISEASON)+RESICE
331 !       IF ((wetflag).OR.(RAINFLAG)) THEN
332 !         IF (ILAND.EQ.1) THEN
333 !           RLUX=50.
334 !         ELSE
335 !           RLUX=100.
336 !         END IF
337 !       END IF
338 !       RCX(P_SO2)=1./(1./RSMX+1./RLUX+1./(RDC+RCLX)
339 !    &                +1./(RAC(ILAND,ISEASON)+RGSX))
340 !       IF (RCX(P_SO2).LT.1.) RCX(P_SO2)=1.
341 
342 !     SO2 according to Erisman et al. 1994
343 !       R_STOM
344         rsmx = rs*dratio(p_so2)
345 !       R_EXT
346         IF (tc>(-1.)) THEN
347           IF (rh<81.3) THEN
348             rlux = 25000.*exp(-0.0693*rh)
349           ELSE
350             rlux = 0.58E12*exp(-0.278*rh)
351           END IF
352         END IF
353         IF (((wetflag) .OR. (rainflag)) .AND. (tc>(-1.))) THEN
354           rlux = 1.
355         END IF
356         IF ((tc>=(-5.)) .AND. (tc<=(-1.))) THEN
357           rlux = 200.
358         END IF
359         IF (tc<(-5.)) THEN
360           rlux = 500.
361         END IF
362 !       INSTEAD OF R_INC R_CL and R_DC of Wesely are used
363         rclx = rcls(iland,iseason)
364 !       DRY SURFACE
365         rgsx = 1000.
366 !       WET SURFACE
367         IF ((wetflag) .OR. (rainflag)) THEN
368           IF (highnh3) THEN
369             rgsx = 0.
370           ELSE
371             rgsx = 500.
372           END IF
373         END IF
374 !       WATER
375         IF (iland==iswater_temp) THEN
376           rgsx = 0.
377         END IF
378 !       SNOW
379         IF ((iseason==4) .OR. (iland==isice_temp)) THEN
380           IF (tc>2.) THEN
381             rgsx = 0.
382           END IF
383           IF ((tc>=(-1.)) .AND. (tc<=2.)) THEN
384             rgsx = 70.*(2.-tc)
385           END IF
386           IF (tc<(-1.)) THEN
387             rgsx = 500.
388           END IF
389         END IF
390 !       TOTAL SURFACE RESISTENCE
391         IF ((iseason/=4) .AND. (ixxxlu(iland)/=1) .AND. (iland/=iswater_temp) .AND. &
392             (iland/=isice_temp)) THEN
393           rcx(p_so2) = 1./(1./rsmx+1./rlux+1./(rclx+rdc+rgsx))
394         ELSE
395           rcx(p_so2) = rgsx
396         END IF
397         IF (rcx(p_so2)<1.) rcx(p_so2) = 1.
398 !     NH3 according to Erisman et al. 1994
399 !       R_STOM
400         rsmx = rs*dratio(p_nh3)
401 !       GRASSLAND (PASTURE DURING GRAZING)
402         IF (ixxxlu(iland)==3) THEN
403           IF (iseason==1) THEN
404 !           SUMMER
405             rcx(p_nh3) = 1000.
406           END IF
407           IF ((iseason==2) .OR. (iseason==3) .OR. (iseason==5)) THEN
408 !           WINTER, NO SNOW
409             IF (tc>-1.) THEN
410               IF (rad/=0.) THEN
411                 rcx(p_nh3) = 50.
412               ELSE
413                 rcx(p_nh3) = 100.
414               END IF
415               IF ((wetflag) .OR. (rainflag)) THEN
416                 rcx(p_nh3) = 20.
417               END IF
418             END IF
419             IF ((tc>=(-5.)) .AND. (tc<=-1.)) THEN
420               rcx(p_nh3) = 200.
421             END IF
422             IF (tc<(-5.)) THEN
423               rcx(p_nh3) = 500.
424             END IF
425           END IF
426         END IF
427 !       AGRICULTURAL LAND (CROPS AND UNGRAZED PASTURE)
428         IF (ixxxlu(iland)==2) THEN
429           IF (iseason==1) THEN
430 !           SUMMER
431             IF (rad/=0.) THEN
432               rcx(p_nh3) = rsmx
433             ELSE
434               rcx(p_nh3) = 200.
435             END IF
436             IF ((wetflag) .OR. (rainflag)) THEN
437               rcx(p_nh3) = 50.
438             END IF
439           END IF
440           IF ((iseason==2) .OR. (iseason==3) .OR. (iseason==5)) THEN
441 !           WINTER, NO SNOW
442             IF (tc>-1.) THEN
443               IF (rad/=0.) THEN
444                 rcx(p_nh3) = rsmx
445               ELSE
446                 rcx(p_nh3) = 300.
447               END IF
448               IF ((wetflag) .OR. (rainflag)) THEN
449                 rcx(p_nh3) = 100.
450               END IF
451             END IF
452             IF ((tc>=(-5.)) .AND. (tc<=-1.)) THEN
453               rcx(p_nh3) = 200.
454             END IF
455             IF (tc<(-5.)) THEN
456               rcx(p_nh3) = 500.
457             END IF
458           END IF
459         END IF
460 !       SEMI-NATURAL ECOSYSTEMS AND FORESTS
461         IF ((ixxxlu(iland)==4) .OR. (ixxxlu(iland)==5) .OR. (ixxxlu( &
462             iland)==6)) THEN
463           IF (rad/=0.) THEN
464             rcx(p_nh3) = 500.
465           ELSE
466             rcx(p_nh3) = 1000.
467           END IF
468           IF ((wetflag) .OR. (rainflag)) THEN
469             IF (highnh3) THEN
470               rcx(p_nh3) = 100.
471             ELSE
472               rcx(p_nh3) = 0.
473             END IF
474           END IF
475           IF ((iseason==2) .OR. (iseason==3) .OR. (iseason==5)) THEN
476 !           WINTER, NO SNOW
477             IF ((tc>=(-5.)) .AND. (tc<=-1.)) THEN
478               rcx(p_nh3) = 200.
479             END IF
480             IF (tc<(-5.)) THEN
481               rcx(p_nh3) = 500.
482             END IF
483           END IF
484         END IF
485 !       WATER
486         IF (iland==iswater_temp) THEN
487           rcx(p_nh3) = 0.
488         END IF
489 !       URBAN AND DESERT (SOIL SURFACES)
490         IF (ixxxlu(iland)==1) THEN
491           IF ( .NOT. wetflag) THEN
492             rcx(p_nh3) = 50.
493           ELSE
494             rcx(p_nh3) = 0.
495           END IF
496         END IF
497 !       SNOW COVERED SURFACES OR PERMANENT ICE
498         IF ((iseason==4) .OR. (iland==isice_temp)) THEN
499           IF (tc>2.) THEN
500             rcx(p_nh3) = 0.
501           END IF
502           IF ((tc>=(-1.)) .AND. (tc<=2.)) THEN
503             rcx(p_nh3) = 70.*(2.-tc)
504           END IF
505           IF (tc<(-1.)) THEN
506             rcx(p_nh3) = 500.
507           END IF
508         END IF
509         IF (rcx(p_nh3)<1.) rcx(p_nh3) = 1.
510 
511       END SUBROUTINE rc
512 ! **********************************************************************
513 ! ************************ SUBROUTINE DEPPART **************************
514 ! **********************************************************************
515       SUBROUTINE deppart(rmol,ustar,rh,clw,iland,dvpart,dvfog)
516 !     THIS SUBROUTINE CALCULATES SURFACE DEPOSITION VELOCITIES
517 !     FOR FINE AEROSOL PARTICLES ACCORDING TO THE MODEL OF
518 !     J. W. ERISMAN, A. VAN PUL, AND P. WYERS,
519 !     ATMOSPHERIC ENVIRONMENT 28 (1994), 2595-2607
520 !     WRITTEN BY WINFRIED SEIDL, APRIL 1997
521 !     MODIFIED BY WINFRIED SEIDL, MARCH 2000
522 !            FOR MM5 VERSION 3
523 ! ----------------------------------------------------------------------
524 ! .. Scalar Arguments ..
525         REAL :: clw, dvfog, dvpart, rh, rmol, ustar
526         INTEGER :: iland
527 ! ..
528 ! .. Intrinsic Functions ..
529         INTRINSIC exp
530 ! ..
531         dvpart = ustar/kpart(iland)
532         IF (rmol<0.) THEN
533 !         INSTABLE LAYERING CORRECTION
534           dvpart = dvpart*(1.+(-300.*rmol)**0.66667)
535         END IF
536         IF (rh>80.) THEN
537 !         HIGH RELATIVE HUMIDITY CORRECTION
538 !         ACCORDING TO J. W. ERISMAN ET AL.
539 !         ATMOSPHERIC ENVIRONMENT 31 (1997), 321-332
540           dvpart = dvpart*(1.+0.37*exp((rh-80.)/20.))
541         END IF
542 
543 !       SEDIMENTATION VELOCITY OF FOG WATER ACCORDING TO
544 !       R. FORKEL, W. SEIDL, R. DLUGI AND E. DEIGELE
545 !       J. GEOPHYS. RES. 95D (1990), 18501-18515
546         dvfog = 0.06*clw
547         IF (ixxxlu(iland)==5) THEN
548 !         TURBULENT DEPOSITION OF FOG WATER IN CONIFEROUS FOREST ACCORDI
549 !         A. T. VERMEULEN ET AL.
550 !         ATMOSPHERIC ENVIRONMENT 31 (1997), 375-386
551           dvfog = dvfog + 0.195*ustar*ustar
552         END IF
553 
554       END SUBROUTINE deppart
555       SUBROUTINE landusevg(vgs,ustar,rmol,z0,zz,dvparx,iland,numgas, &
556           srfres,aer_res,p_sulf)
557 !     This subroutine calculates the species specific deposition velocit
558 !     as a function of the local meteorology and land use.  The depositi
559 !     Velocity is also landuse specific.
560 !     Reference: Hsieh, C.M., Wesely, M.L. and Walcek, C.J. (1986)
561 !                A Dry Deposition Module for Regional Acid Deposition
562 !                EPA report under agreement DW89930060-01
563 !     Revised version by Darrell Winner (January 1991)
564 !        Environmental Engineering Science 138-78
565 !           California Institute of Technology
566 !              Pasadena, CA  91125
567 !     Modified by Winfried Seidl (August 1997)
568 !       Fraunhofer-Institut fuer Atmosphaerische Umweltforschung
569 !                    Garmisch-Partenkirchen, D-82467
570 !          for use of Wesely and Erisman surface resistances
571 !     Inputs:
572 !        Ustar  : The grid average friction velocity (m/s)
573 !        Rmol   : Reciprocal of the Monin-Obukhov length (1/m)
574 !        Z0     : Surface roughness height for the grid square (m)
575 !        SrfRes : Array of landuse/atmospheric/species resistances (s/m)
576 !        Slist  : Array of chemical species codes
577 !        Dvparx : Array of surface deposition velocity of fine aerosol p
578 !     Outputs:
579 !        Vgs    : Array of species and landuse specific deposition
580 !                 velocities (m/s)
581 !        Vg     : Cell-average deposition velocity by species (m/s)
582 !     Variables used:
583 !        SCPR23  : (Schmidt #/Prandtl #)**(2/3) Diffusion correction fac
584 !        Zr      : Reference Height (m)
585 !        Iatmo   : Parameter specifying the stabilty class (Function of
586 !        Z0      : Surface roughness height (m)
587 !        karman  : Von Karman constant (from module_model_constants)
588         USE module_model_constants, only: karman
589 !     Local Variables
590 ! .. Scalar Arguments ..
591         REAL :: dvparx, rmol, ustar, z0, zz
592         real :: aer_res, polint
593         INTEGER :: iland, numgas, p_sulf
594 ! ..
595 ! .. Array Arguments ..
596         REAL :: srfres(numgas), vgs(numgas)
597 ! ..
598 ! .. Local Scalars ..
599         REAL :: vgp, vgpart, zr
600         INTEGER :: jspec
601 ! ..
602 ! .. Local Arrays ..
603         REAL :: vgspec(numgas)
604 ! ..
605 !   Set the reference height (10.0 m)
606         zr = 10.0
607 
608 !   CALCULATE THE DEPOSITION VELOCITY without any surface
609 !   resistance term, i.e. 1 / (ra + rb)
610         CALL depvel(numgas,rmol,zr,z0,ustar,vgspec,vgpart,aer_res)
611 
612 !   Calculate the deposition velocity for each species
613 !   and grid cell by looping through all the possibile combinations
614 !   of the two
615 
616         vgp = 1.0/((1.0/vgpart)+(1.0/dvparx))
617 
618 !   Loop through the various species
619 
620         DO jspec = 1, numgas
621 
622 !   Add in the surface resistance term, rc (SrfRes)
623 
624           vgs(jspec) = 1.0/(1.0/vgspec(jspec)+srfres(jspec))
625         END DO
626         vgs(p_sulf) = vgp
627 
628         CALL cellvg(vgs,ustar,zz,zr,rmol,numgas)
629 
630         RETURN
631       END SUBROUTINE landusevg
632 
633       SUBROUTINE cellvg(vgtemp,ustar,dz,zr,rmol,nspec)
634 !     THIS PROGRAM HAS BEEN DESIGNED TO CALCULATE THE CELL AVERAGE
635 !     DEPOSITION VELOCITY GIVEN THE VALUE OF VG AT SOME REFERENCE
636 !     HEIGHT ZR WHICH IS MUCH SMALLER THAN THE CELL HEIGHT DZ.
637 !       PROGRAM WRITTEN BY GREGORY J.MCRAE (NOVEMBER 1977)
638 !         Modified by Darrell A. Winner    (February 1991)
639 !.....PROGRAM VARIABLES...
640 !     VgTemp   - DEPOSITION VELOCITY AT THE REFERENCE HEIGHT
641 !     USTAR    - FRICTION VELOCITY
642 !     RMOL     - RECIPROCAL OF THE MONIN-OBUKHOV LENGTH
643 !     ZR       - REFERENCE HEIGHT
644 !     DZ       - CELL HEIGHT
645 !     CELLVG   - CELL AVERAGE DEPOSITION VELOCITY
646 !     VK       - VON KARMAN CONSTANT
647         USE module_model_constants, only: karman
648 !     Local Variables
649 ! .. Scalar Arguments ..
650         REAL :: dz, rmol, ustar, zr
651         INTEGER :: nspec
652 ! ..
653 ! .. Array Arguments ..
654         REAL :: vgtemp(nspec)
655 ! ..
656 ! .. Local Scalars ..
657         REAL :: a, fac, pdz, pzr, vk
658         INTEGER :: nss
659 ! ..
660 ! .. Intrinsic Functions ..
661         INTRINSIC alog, sqrt
662 ! ..
663 !     Set the von Karman constant
664         vk = karman
665 
666 !     DETERMINE THE STABILITY BASED ON THE CONDITIONS
667 !             1/L < 0 UNSTABLE
668 !             1/L = 0 NEUTRAL
669 !             1/L > 0 STABLE
670 
671 
672         DO nss = 1, nspec
673           IF (rmol<0) THEN
674             pdz = sqrt(1.0-9.0*dz*rmol)
675             pzr = sqrt(1.0-9.0*zr*rmol)
676             fac = ((pdz-1.0)/(pzr-1.0))*((pzr+1.0)/(pdz+1.0))
677             a = 0.74*dz*alog(fac) + (0.164/rmol)*(pdz-pzr)
678           ELSE IF (rmol==0) THEN
679             a = 0.74*(dz*alog(dz/zr)-dz+zr)
680           ELSE
681             a = 0.74*(dz*alog(dz/zr)-dz+zr) + (2.35*rmol)*(dz-zr)**2
682           END IF
683 
684 !     CALCULATE THE DEPOSITION VELOCITIY
685 
686           vgtemp(nss) = vgtemp(nss)/(1.0+vgtemp(nss)*a/(vk*ustar*(dz-zr)))
687         END DO
688 
689         RETURN
690       END SUBROUTINE cellvg
691       SUBROUTINE depvel(numgas,rmol,zr,z0,ustar,depv,vgpart,aer_res)
692 !     THIS FUNCTION HAS BEEN DESIGNED TO EVALUATE AN UPPER LIMIT
693 !     FOR THE POLLUTANT DEPOSITION VELOCITY AS A FUNCTION OF THE
694 !     SURFACE ROUGHNESS AND METEOROLOGICAL CONDITIONS.
695 !     PROGRAM WRITTEN BY GREGORY J.MCRAE (NOVEMBER 1977)
696 !         Modified by Darrell A. Winner  (Feb. 1991)
697 !                  by Winfried Seidl     (Aug. 1997)
698 !.....PROGRAM VARIABLES...
699 !     RMOL     - RECIPROCAL OF THE MONIN-OBUKHOV LENGTH
700 !     ZR       - REFERENCE HEIGHT
701 !     Z0       - SURFACE ROUGHNESS HEIGHT
702 !     SCPR23   - (Schmidt #/Prandtl #)**(2/3) Diffusion correction fact
703 !     UBAR     - ABSOLUTE VALUE OF SURFACE WIND SPEED
704 !     DEPVEL   - POLLUTANT DEPOSITION VELOCITY
705 !     Vk       - VON KARMAN CONSTANT
706 !     USTAR    - FRICTION VELOCITY U*
707 !     POLINT   - POLLUTANT INTEGRAL
708 !     AER_RES  - AERODYNAMIC RESISTANCE
709 !.....REFERENCES...
710 !     MCRAE, G.J. ET AL. (1983) MATHEMATICAL MODELING OF PHOTOCHEMICAL
711 !       AIR POLLUTION, ENVIRONMENTAL QUALITY LABORATORY REPORT 18,
712 !       CALIFORNIA INSTITUTE OF TECHNOLOGY, PASADENA, CALIFORNIA.
713 !.....RESTRICTIONS...
714 !     1. THE MODEL EDDY DIFFUSIVITIES ARE BASED ON MONIN-OBUKHOV
715 !        SIMILARITY THEORY AND SO ARE ONLY APPLICABLE IN THE
716 !        SURFACE LAYER, A HEIGHT OF O(30M).
717 !     2. ALL INPUT UNITS MUST BE CONSISTENT
718 !     3. THE PHI FUNCTIONS USED TO CALCULATE THE FRICTION
719 !        VELOCITY U* AND THE POLLUTANT INTEGRALS ARE BASED
720 !        ON THE WORK OF BUSINGER ET AL.(1971).
721 !     4. THE MOMENTUM AND POLLUTANT DIFFUSIVITIES ARE NOT
722 !        THE SAME FOR THE CASES L<0 AND L>0.
723         USE module_model_constants, only: karman
724 !     Local Variables
725 ! .. Scalar Arguments ..
726         REAL :: rmol, ustar, vgpart, z0, zr, aer_res
727         INTEGER :: numgas
728 ! ..
729 ! .. Array Arguments ..
730         REAL :: depv(numgas)
731 ! ..
732 ! .. Local Scalars ..
733         REAL :: ao, ar, polint, vk
734         INTEGER :: l
735 ! ..
736 ! .. Intrinsic Functions ..
737         INTRINSIC alog
738 ! ..
739 !     Set the von Karman constant
740         vk = karman
741 
742 !     Calculate the diffusion correction factor
743 !     SCPR23 is calculated as (Sc/Pr)**(2/3) using Sc= 1.15 and Pr= 1.0
744 !     SCPR23 = 1.10
745 
746 !     DETERMINE THE STABILITY BASED ON THE CONDITIONS
747 !             1/L < 0 UNSTABLE
748 !             1/L = 0 NEUTRAL
749 !             1/L > 0 STABLE
750 
751         if(abs(rmol) < 1.E-6 ) rmol = 0.
752 
753         IF (rmol<0) THEN
754           ar = ((1.0-9.0*zr*rmol)**(0.25)+0.001)**2
755           ao = ((1.0-9.0*z0*rmol)**(0.25)+0.001)**2
756           polint = 0.74*(alog((ar-1.0)/(ar+1.0))-alog((ao-1.0)/(ao+1.0)))
757         ELSE IF (rmol==0.) THEN
758           polint = 0.74*alog(zr/z0)
759         ELSE
760           polint = 0.74*alog(zr/z0) + 4.7*rmol*(zr-z0)
761         END IF
762 
763 !     CALCULATE THE Maximum DEPOSITION VELOCITY
764 
765         DO l = 1, numgas
766           depv(l) = ustar*vk/(2.0*scpr23(l)+polint)
767         END DO
768         vgpart = ustar*vk/polint
769         aer_res = polint/(karman*max(ustar,1.0e-4))
770 
771         RETURN
772       END SUBROUTINE depvel
773 
774       SUBROUTINE dep_init(id,config_flags,numgas)
775   USE module_model_constants
776   USE module_configure
777   USE module_state_description                       
778    TYPE (grid_config_rec_type) , INTENT (in) ::     config_flags
779 ! ..
780 ! .. Scalar Arguments ..
781         integer, intent(in) :: id, numgas
782 
783 ! ..
784 ! .. Local Scalars ..
785         REAL :: sc
786         INTEGER :: iland, iseason, l
787         integer :: iprt
788 ! ..
789 ! .. Local Arrays ..
790         REAL :: dat1(nlu,dep_seasons), dat2(nlu,dep_seasons),         &
791                 dat3(nlu,dep_seasons), dat4(nlu,dep_seasons),         &
792                 dat5(nlu,dep_seasons), dat6(nlu,dep_seasons),         &
793                 dat7(nlu,dep_seasons), dvj(numgas)
794 ! ..
795 ! .. Make sure that the model is being run with a soil model. Otherwise,
796 !    iland will be zero in deppart, which will try to pull non-exisant
797 !    array locations.
798         call nl_get_sf_surface_physics(id,l)
799         if( l == 0 ) &
800              call wrf_error_fatal("ERROR: Cannot use dry deposition without using a soil model.")
801 
802 ! ..
803 ! .. Data Statements ..
804 !     RI for stomatal resistance
805 !      data ((ri(ILAND,ISEASON),ILAND=1,nlu),ISEASON=1,dep_seasons)/0.10E+11, &
806         DATA ((dat1(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+11, &
807           0.60E+02, 0.60E+02, 0.60E+02, 0.60E+02, 0.70E+02, 0.12E+03, &
808           0.12E+03, 0.12E+03, 0.12E+03, 0.70E+02, 0.13E+03, 0.70E+02, &
809           0.13E+03, 0.10E+03, 0.10E+11, 0.80E+02, 0.10E+03, 0.10E+11, &
810           0.80E+02, 0.10E+03, 0.10E+03, 0.10E+11, 0.10E+11, 0.10E+11, &
811           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
812           0.10E+11, 0.10E+11, 0.10E+11, 0.12E+03, 0.10E+11, 0.10E+11, &
813           0.70E+02, 0.25E+03, 0.50E+03, 0.10E+11, 0.10E+11, 0.50E+03, &
814           0.10E+11, 0.10E+11, 0.50E+03, 0.50E+03, 0.10E+11, 0.10E+11, &
815           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
816           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.12E+03, 0.10E+11, &
817           0.10E+11, 0.70E+02, 0.25E+03, 0.50E+03, 0.10E+11, 0.10E+11, &
818           0.50E+03, 0.10E+11, 0.10E+11, 0.50E+03, 0.50E+03, 0.10E+11, &
819           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
820           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
821           0.10E+11, 0.10E+11, 0.70E+02, 0.40E+03, 0.80E+03, 0.10E+11, &
822           0.10E+11, 0.80E+03, 0.10E+11, 0.10E+11, 0.80E+03, 0.80E+03, &
823           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.12E+03, 0.12E+03, &
824           0.12E+03, 0.12E+03, 0.14E+03, 0.24E+03, 0.24E+03, 0.24E+03, &
825           0.12E+03, 0.14E+03, 0.25E+03, 0.70E+02, 0.25E+03, 0.19E+03, &
826           0.10E+11, 0.16E+03, 0.19E+03, 0.10E+11, 0.16E+03, 0.19E+03, &
827           0.19E+03, 0.10E+11, 0.10E+11, 0.10E+11/
828 ! ..
829         IF (nlu/=25) THEN
830           call wrf_debug(0, 'number of land use classifications not correct ')
831           CALL wrf_error_fatal ( "LAND USE CLASSIFICATIONS NOT 25")
832         END IF
833         IF (dep_seasons/=5) THEN
834           call wrf_debug(0, 'number of dep_seasons not correct ')
835           CALL wrf_error_fatal ( "DEP_SEASONS NOT 5")
836         END IF
837 
838 !     SURFACE RESISTANCE DATA FOR DEPOSITION MODEL OF
839 !     M. L. WESELY, ATMOSPHERIC ENVIRONMENT 23 (1989) 1293-1304
840 
841 !     Seasonal categories:
842 !     1: midsummer with lush vegetation
843 !     2: autumn with unharvested cropland
844 !     3: late autumn with frost, no snow
845 !     4: winter, snow on ground and subfreezing
846 !     5: transitional spring with partially green short annuals
847 
848 !     Land use types:
849 !     USGS type                                Wesely type
850 !      1: Urban and built-up land              1
851 !      2: Dryland cropland and pasture         2
852 !      3: Irrigated cropland and pasture       2
853 !      4: Mix. dry/irrg. cropland and pasture  2
854 !      5: Cropland/grassland mosaic            2
855 !      6: Cropland/woodland mosaic             4
856 !      7: Grassland                            3
857 !      8: Shrubland                            3
858 !      9: Mixed shrubland/grassland            3
859 !     10: Savanna                              3, always summer
860 !     11: Deciduous broadleaf forest           4
861 !     12: Deciduous needleleaf forest          5, autumn and winter modi
862 !     13: Evergreen broadleaf forest           4, always summer
863 !     14: Evergreen needleleaf forest          5
864 !     15: Mixed Forest                         6
865 !     16: Water Bodies                         7
866 !     17: Herbaceous wetland                   9
867 !     18: Wooded wetland                       6
868 !     19: Barren or sparsely vegetated         8
869 !     20: Herbaceous Tundra                    9
870 !     21: Wooded Tundra                        6
871 !     22: Mixed Tundra                         6
872 !     23: Bare Ground Tundra                   8
873 !     24: Snow or Ice                          -, always winter
874 !     25: No data                              8
875 
876 
877 !     Order of data:
878 !      |
879 !      |   seasonal category
880 !     \|/
881 !     ---> landuse type
882 !     1       2       3       4       5       6       7       8       9
883 !     RLU for outer surfaces in the upper canopy
884         DO iseason = 1, dep_seasons
885           DO iland = 1, nlu
886             ri(iland,iseason) = dat1(iland,iseason)
887           END DO
888         END DO
889 !      data ((rlu(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.10E+11, &
890         DATA ((dat2(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+11, &
891           0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, &
892           0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, &
893           0.20E+04, 0.20E+04, 0.10E+11, 0.25E+04, 0.20E+04, 0.10E+11, &
894           0.25E+04, 0.20E+04, 0.20E+04, 0.10E+11, 0.10E+11, 0.10E+11, &
895           0.10E+11, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, &
896           0.90E+04, 0.90E+04, 0.90E+04, 0.20E+04, 0.90E+04, 0.90E+04, &
897           0.20E+04, 0.40E+04, 0.80E+04, 0.10E+11, 0.90E+04, 0.80E+04, &
898           0.10E+11, 0.90E+04, 0.80E+04, 0.80E+04, 0.10E+11, 0.10E+11, &
899           0.10E+11, 0.10E+11, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, &
900           0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, 0.20E+04, 0.90E+04, &
901           0.90E+04, 0.20E+04, 0.40E+04, 0.80E+04, 0.10E+11, 0.90E+04, &
902           0.80E+04, 0.10E+11, 0.90E+04, 0.80E+04, 0.80E+04, 0.10E+11, &
903           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
904           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
905           0.10E+11, 0.10E+11, 0.20E+04, 0.60E+04, 0.90E+04, 0.10E+11, &
906           0.90E+04, 0.90E+04, 0.10E+11, 0.90E+04, 0.90E+04, 0.90E+04, &
907           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.40E+04, 0.40E+04, &
908           0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, &
909           0.20E+04, 0.40E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.30E+04, &
910           0.10E+11, 0.40E+04, 0.30E+04, 0.10E+11, 0.40E+04, 0.30E+04, &
911           0.30E+04, 0.10E+11, 0.10E+11, 0.10E+11/
912         DO iseason = 1, dep_seasons
913           DO iland = 1, nlu
914             rlu(iland,iseason) = dat2(iland,iseason)
915           END DO
916         END DO
917 !     RAC for transfer that depends on canopy height and density
918 !      data ((rac(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.10E+03, &
919         DATA ((dat3(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+03, &
920           0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+04, 0.10E+03, &
921           0.10E+03, 0.10E+03, 0.10E+03, 0.20E+04, 0.20E+04, 0.20E+04, &
922           0.20E+04, 0.20E+04, 0.00E+00, 0.30E+03, 0.20E+04, 0.00E+00, &
923           0.30E+03, 0.20E+04, 0.20E+04, 0.00E+00, 0.00E+00, 0.00E+00, &
924           0.10E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+04, &
925           0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, 0.15E+04, 0.20E+04, &
926           0.20E+04, 0.20E+04, 0.17E+04, 0.00E+00, 0.20E+03, 0.17E+04, &
927           0.00E+00, 0.20E+03, 0.17E+04, 0.17E+04, 0.00E+00, 0.00E+00, &
928           0.00E+00, 0.10E+03, 0.10E+02, 0.10E+02, 0.10E+02, 0.10E+02, &
929           0.10E+04, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+04, &
930           0.20E+04, 0.20E+04, 0.20E+04, 0.15E+04, 0.00E+00, 0.10E+03, &
931           0.15E+04, 0.00E+00, 0.10E+03, 0.15E+04, 0.15E+04, 0.00E+00, &
932           0.00E+00, 0.00E+00, 0.10E+03, 0.10E+02, 0.10E+02, 0.10E+02, &
933           0.10E+02, 0.10E+04, 0.10E+02, 0.10E+02, 0.10E+02, 0.10E+02, &
934           0.10E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.15E+04, 0.00E+00, &
935           0.50E+02, 0.15E+04, 0.00E+00, 0.50E+02, 0.15E+04, 0.15E+04, &
936           0.00E+00, 0.00E+00, 0.00E+00, 0.10E+03, 0.50E+02, 0.50E+02, &
937           0.50E+02, 0.50E+02, 0.12E+04, 0.80E+02, 0.80E+02, 0.80E+02, &
938           0.10E+03, 0.12E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.15E+04, &
939           0.00E+00, 0.20E+03, 0.15E+04, 0.00E+00, 0.20E+03, 0.15E+04, &
940           0.15E+04, 0.00E+00, 0.00E+00, 0.00E+00/
941         DO iseason = 1, dep_seasons
942           DO iland = 1, nlu
943             rac(iland,iseason) = dat3(iland,iseason)
944           END DO
945         END DO
946 !     RGSS for ground surface  SO2
947 !      data ((rgss(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.40E+03, &
948         DATA ((dat4(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.40E+03, &
949           0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.50E+03, 0.35E+03, &
950           0.35E+03, 0.35E+03, 0.35E+03, 0.50E+03, 0.50E+03, 0.50E+03, &
951           0.50E+03, 0.10E+03, 0.10E+01, 0.10E+01, 0.10E+03, 0.10E+04, &
952           0.10E+01, 0.10E+03, 0.10E+03, 0.10E+04, 0.10E+03, 0.10E+04, &
953           0.40E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.50E+03, &
954           0.35E+03, 0.35E+03, 0.35E+03, 0.35E+03, 0.50E+03, 0.50E+03, &
955           0.50E+03, 0.50E+03, 0.10E+03, 0.10E+01, 0.10E+01, 0.10E+03, &
956           0.10E+04, 0.10E+01, 0.10E+03, 0.10E+03, 0.10E+04, 0.10E+03, &
957           0.10E+04, 0.40E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, &
958           0.50E+03, 0.35E+03, 0.35E+03, 0.35E+03, 0.35E+03, 0.50E+03, &
959           0.50E+03, 0.50E+03, 0.50E+03, 0.20E+03, 0.10E+01, 0.10E+01, &
960           0.20E+03, 0.10E+04, 0.10E+01, 0.20E+03, 0.20E+03, 0.10E+04, &
961           0.10E+03, 0.10E+04, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, &
962           0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, 0.10E+03, &
963           0.10E+03, 0.10E+03, 0.50E+03, 0.10E+03, 0.10E+03, 0.10E+01, &
964           0.10E+03, 0.10E+03, 0.10E+04, 0.10E+03, 0.10E+03, 0.10E+03, &
965           0.10E+04, 0.10E+03, 0.10E+04, 0.50E+03, 0.15E+03, 0.15E+03, &
966           0.15E+03, 0.15E+03, 0.50E+03, 0.35E+03, 0.35E+03, 0.35E+03, &
967           0.35E+03, 0.50E+03, 0.50E+03, 0.50E+03, 0.50E+03, 0.20E+03, &
968           0.10E+01, 0.10E+01, 0.20E+03, 0.10E+04, 0.10E+01, 0.20E+03, &
969           0.20E+03, 0.10E+04, 0.10E+03, 0.10E+04/
970         DO iseason = 1, dep_seasons
971           DO iland = 1, nlu
972             rgss(iland,iseason) = dat4(iland,iseason)
973           END DO
974         END DO
975 !     RGSO for ground surface  O3
976 !      data ((rgso(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.30E+03, &
977         DATA ((dat5(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.30E+03, &
978           0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.20E+03, 0.20E+03, &
979           0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, &
980           0.20E+03, 0.30E+03, 0.20E+04, 0.10E+04, 0.30E+03, 0.40E+03, &
981           0.10E+04, 0.30E+03, 0.30E+03, 0.40E+03, 0.35E+04, 0.40E+03, &
982           0.30E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.20E+03, &
983           0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, &
984           0.20E+03, 0.20E+03, 0.30E+03, 0.20E+04, 0.80E+03, 0.30E+03, &
985           0.40E+03, 0.80E+03, 0.30E+03, 0.30E+03, 0.40E+03, 0.35E+04, &
986           0.40E+03, 0.30E+03, 0.15E+03, 0.15E+03, 0.15E+03, 0.15E+03, &
987           0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, &
988           0.20E+03, 0.20E+03, 0.20E+03, 0.30E+03, 0.20E+04, 0.10E+04, &
989           0.30E+03, 0.40E+03, 0.10E+04, 0.30E+03, 0.30E+03, 0.40E+03, &
990           0.35E+04, 0.40E+03, 0.60E+03, 0.35E+04, 0.35E+04, 0.35E+04, &
991           0.35E+04, 0.35E+04, 0.35E+04, 0.35E+04, 0.35E+04, 0.35E+04, &
992           0.35E+04, 0.35E+04, 0.20E+03, 0.35E+04, 0.35E+04, 0.20E+04, &
993           0.35E+04, 0.35E+04, 0.40E+03, 0.35E+04, 0.35E+04, 0.35E+04, &
994           0.40E+03, 0.35E+04, 0.40E+03, 0.30E+03, 0.15E+03, 0.15E+03, &
995           0.15E+03, 0.15E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, &
996           0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.20E+03, 0.30E+03, &
997           0.20E+04, 0.10E+04, 0.30E+03, 0.40E+03, 0.10E+04, 0.30E+03, &
998           0.30E+03, 0.40E+03, 0.35E+04, 0.40E+03/
999         DO iseason = 1, dep_seasons
1000           DO iland = 1, nlu
1001             rgso(iland,iseason) = dat5(iland,iseason)
1002           END DO
1003         END DO
1004 !     RCLS for exposed surfaces in the lower canopy  SO2
1005 !      data ((rcls(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.10E+11, &
1006         DATA ((dat6(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+11, &
1007           0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, &
1008           0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.20E+04, &
1009           0.20E+04, 0.20E+04, 0.10E+11, 0.25E+04, 0.20E+04, 0.10E+11, &
1010           0.25E+04, 0.20E+04, 0.20E+04, 0.10E+11, 0.10E+11, 0.10E+11, &
1011           0.10E+11, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, &
1012           0.90E+04, 0.90E+04, 0.90E+04, 0.20E+04, 0.90E+04, 0.90E+04, &
1013           0.20E+04, 0.20E+04, 0.40E+04, 0.10E+11, 0.90E+04, 0.40E+04, &
1014           0.10E+11, 0.90E+04, 0.40E+04, 0.40E+04, 0.10E+11, 0.10E+11, &
1015           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
1016           0.90E+04, 0.90E+04, 0.90E+04, 0.90E+04, 0.20E+04, 0.90E+04, &
1017           0.90E+04, 0.20E+04, 0.30E+04, 0.60E+04, 0.10E+11, 0.90E+04, &
1018           0.60E+04, 0.10E+11, 0.90E+04, 0.60E+04, 0.60E+04, 0.10E+11, &
1019           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
1020           0.10E+11, 0.90E+04, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, &
1021           0.90E+04, 0.90E+04, 0.20E+04, 0.20E+03, 0.40E+03, 0.10E+11, &
1022           0.90E+04, 0.40E+03, 0.10E+11, 0.90E+04, 0.40E+03, 0.40E+03, &
1023           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.40E+04, 0.40E+04, &
1024           0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, 0.40E+04, &
1025           0.20E+04, 0.40E+04, 0.20E+04, 0.20E+04, 0.20E+04, 0.30E+04, &
1026           0.10E+11, 0.40E+04, 0.30E+04, 0.10E+11, 0.40E+04, 0.30E+04, &
1027           0.30E+04, 0.10E+11, 0.10E+11, 0.10E+11/
1028         DO iseason = 1, dep_seasons
1029           DO iland = 1, nlu
1030             rcls(iland,iseason) = dat6(iland,iseason)
1031           END DO
1032         END DO
1033 !     RCLO for exposed surfaces in the lower canopy  O3
1034 !      data ((rclo(ILAND,ISEASON),ILAND=1,25),ISEASON=1,5)/0.10E+11, &
1035         DATA ((dat7(iland,iseason),iland=1,nlu),iseason=1,dep_seasons)/0.10E+11, &
1036           0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, &
1037           0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, &
1038           0.10E+04, 0.10E+04, 0.10E+11, 0.10E+04, 0.10E+04, 0.10E+11, &
1039           0.10E+04, 0.10E+04, 0.10E+04, 0.10E+11, 0.10E+11, 0.10E+11, &
1040           0.10E+11, 0.40E+03, 0.40E+03, 0.40E+03, 0.40E+03, 0.40E+03, &
1041           0.40E+03, 0.40E+03, 0.40E+03, 0.10E+04, 0.40E+03, 0.40E+03, &
1042           0.10E+04, 0.10E+04, 0.60E+03, 0.10E+11, 0.40E+03, 0.60E+03, &
1043           0.10E+11, 0.40E+03, 0.60E+03, 0.60E+03, 0.10E+11, 0.10E+11, &
1044           0.10E+11, 0.10E+11, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, &
1045           0.40E+03, 0.40E+03, 0.40E+03, 0.40E+03, 0.10E+04, 0.40E+03, &
1046           0.40E+03, 0.10E+04, 0.10E+04, 0.60E+03, 0.10E+11, 0.80E+03, &
1047           0.60E+03, 0.10E+11, 0.80E+03, 0.60E+03, 0.60E+03, 0.10E+11, &
1048           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+04, 0.10E+04, 0.10E+04, &
1049           0.10E+04, 0.40E+03, 0.10E+04, 0.10E+04, 0.10E+04, 0.10E+04, &
1050           0.40E+03, 0.40E+03, 0.10E+04, 0.15E+04, 0.60E+03, 0.10E+11, &
1051           0.80E+03, 0.60E+03, 0.10E+11, 0.80E+03, 0.60E+03, 0.60E+03, &
1052           0.10E+11, 0.10E+11, 0.10E+11, 0.10E+11, 0.10E+04, 0.10E+04, &
1053           0.10E+04, 0.10E+04, 0.50E+03, 0.50E+03, 0.50E+03, 0.50E+03, &
1054           0.10E+04, 0.50E+03, 0.15E+04, 0.10E+04, 0.15E+04, 0.70E+03, &
1055           0.10E+11, 0.60E+03, 0.70E+03, 0.10E+11, 0.60E+03, 0.70E+03, &
1056           0.70E+03, 0.10E+11, 0.10E+11, 0.10E+11/
1057         DO iseason = 1, dep_seasons
1058           DO iland = 1, nlu
1059             rclo(iland,iseason) = dat7(iland,iseason)
1060           END DO
1061         END DO
1062         DO l = 1, numgas
1063           hstar(l) = 0.
1064           hstar4(l) = 0.
1065           dhr(l) = 0.
1066           f0(l) = 0.
1067           dvj(l) = 99.
1068         END DO
1069 
1070 !     HENRY''S LAW COEFFICIENTS
1071 !     Effective Henry''s law coefficient at pH 7
1072 !     [KH298]=mole/(l atm)
1073 
1074         hstar(p_no2) = 6.40E-3
1075         hstar(p_no) = 1.90E-3
1076         hstar(p_pan) = 2.97E+0
1077         hstar(p_o3) = 1.13E-2
1078         hstar(p_hcho) = 2.97E+3
1079         hstar(p_aco3) = 1.14E+1
1080         hstar(p_tpan) = 2.97E+0
1081         hstar(p_hono) = 3.47E+5
1082         hstar(p_no3) = 1.50E+1
1083         hstar(p_hno4) = 2.00E+13
1084         hstar(p_h2o2) = 7.45E+4
1085         hstar(p_co) = 8.20E-3
1086         hstar(p_ald) = 1.14E+1
1087         hstar(p_op1) = 2.21E+2
1088         hstar(p_op2) = 1.68E+6
1089         hstar(p_paa) = 4.73E+2
1090         hstar(p_ket) = 3.30E+1
1091         hstar(p_gly) = 1.40E+6
1092         hstar(p_mgly) = 3.71E+3
1093         hstar(p_dcb) = 1.40E+6
1094         hstar(p_onit) = 1.13E+0
1095         hstar(p_so2) = 2.53E+5
1096         hstar(p_eth) = 2.00E-3
1097         hstar(p_hc3) = 1.42E-3
1098         hstar(p_hc5) = 1.13E-3
1099         hstar(p_hc8) = 1.42E-3
1100         hstar(p_olt) = 4.76E-3
1101         hstar(p_oli) = 1.35E-3
1102         hstar(p_tol) = 1.51E-1
1103         hstar(p_csl) = 4.00E+5
1104         hstar(p_xyl) = 1.45E-1
1105         hstar(p_iso) = 4.76E-3
1106         hstar(p_hno3) = 2.69E+13
1107         hstar(p_ora1) = 9.85E+6
1108         hstar(p_ora2) = 9.63E+5
1109         hstar(p_nh3) = 1.04E+4
1110         hstar(p_n2o5) = 1.00E+10
1111         if(p_ol2.gt.1)hstar(p_ol2) = 4.67E-3
1112         hstar(p_par) = 1.13E-3  !wig, 1-May-2007: for CBMZ
1113 !
1114 !  FOLLOWING FOR RACM
1115 !
1116         if(p_ete.gt.1)then
1117            HSTAR(p_ETE )=4.67E-3
1118            HSTAR(p_API )=4.76E-3
1119            HSTAR(p_LIM )=4.76E-3
1120            HSTAR(p_DIEN)=4.76E-3
1121            HSTAR(p_CH4 )=1.50E-3
1122            HSTAR(p_CO2 )=1.86E-1
1123            HSTAR(p_MACR)=1.14E+1
1124            HSTAR(p_UDD )=1.40E+6
1125            HSTAR(p_HKET)=7.80E+3
1126            DHR(p_ETE )=    0.
1127            DHR(p_API )=    0.
1128            DHR(p_LIM )=    0.
1129            DHR(p_DIEN)=    0.
1130            DHR(p_CH4 )=    0.
1131            DHR(p_CO2 )= 1636.
1132            DHR(p_MACR)= 6266.
1133            DHR(p_UDD )=    0.
1134            DHR(p_HKET)=    0.
1135            F0(p_ETE )=0.
1136            F0(p_API )=0.
1137            F0(p_LIM )=0.
1138            F0(p_DIEN)=0.
1139            F0(p_CH4 )=0.
1140            F0(p_CO2 )=0.
1141            F0(p_MACR)=0.
1142            F0(p_UDD )=0.
1143            F0(p_HKET)=0.
1144            DVJ(p_ETE )=0.189
1145            DVJ(p_API )=0.086
1146            DVJ(p_LIM )=0.086
1147            DVJ(p_DIEN)=0.136
1148            DVJ(p_CH4 )=0.250
1149            DVJ(p_CO2 )=0.151
1150            DVJ(p_MACR)=0.120
1151            DVJ(p_UDD )=0.092
1152            DVJ(p_HKET)=0.116
1153         endif
1154 !     -DH/R (for temperature correction)
1155 !     [-DH/R]=K
1156 
1157         dhr(p_no2) = 2500.
1158         dhr(p_no) = 1480.
1159         dhr(p_pan) = 5760.
1160         dhr(p_o3) = 2300.
1161         dhr(p_hcho) = 7190.
1162         dhr(p_aco3) = 6266.
1163         dhr(p_tpan) = 5760.
1164         dhr(p_hono) = 3775.
1165         dhr(p_no3) = 0.
1166         dhr(p_hno4) = 0.
1167         dhr(p_h2o2) = 6615.
1168         dhr(p_co) = 0.
1169         dhr(p_ald) = 6266.
1170         dhr(p_op1) = 5607.
1171         dhr(p_op2) = 10240.
1172         dhr(p_paa) = 6170.
1173         dhr(p_ket) = 5773.
1174         dhr(p_gly) = 0.
1175         dhr(p_mgly) = 7541.
1176         dhr(p_dcb) = 0.
1177         dhr(p_onit) = 5487.
1178         dhr(p_so2) = 5816.
1179         dhr(p_eth) = 0.
1180         dhr(p_hc3) = 0.
1181         dhr(p_hc5) = 0.
1182         dhr(p_hc8) = 0.
1183         dhr(p_olt) = 0.
1184         dhr(p_oli) = 0.
1185         dhr(p_tol) = 0.
1186         dhr(p_csl) = 0.
1187         dhr(p_xyl) = 0.
1188         dhr(p_iso) = 0.
1189         dhr(p_hno3) = 8684.
1190         dhr(p_ora1) = 5716.
1191         dhr(p_ora2) = 8374.
1192         dhr(p_nh3) = 3660.
1193         dhr(p_n2o5) = 0.
1194         if(p_ol2.gt.1)dhr(p_ol2) = 0.
1195         dhr(p_par) = 0.   !wig, 1-May-2007: for CBMZ
1196 !     REACTIVITY FACTORS
1197 !     [f0]=1
1198 
1199         f0(p_no2) = 0.1
1200         f0(p_no) = 0.
1201         f0(p_pan) = 0.1
1202         f0(p_o3) = 1.
1203         f0(p_hcho) = 0.
1204         f0(p_aco3) = 1.
1205         f0(p_tpan) = 0.1
1206         f0(p_hono) = 0.1
1207         f0(p_no3) = 1.
1208         f0(p_hno4) = 0.1
1209         f0(p_h2o2) = 1.
1210         f0(p_co) = 0.
1211         f0(p_ald) = 0.
1212         f0(p_op1) = 0.1
1213         f0(p_op2) = 0.1
1214         f0(p_paa) = 0.1
1215         f0(p_ket) = 0.
1216         f0(p_gly) = 0.
1217         f0(p_mgly) = 0.
1218         f0(p_dcb) = 0.
1219         f0(p_onit) = 0.
1220         f0(p_so2) = 0.
1221         f0(p_eth) = 0.
1222         f0(p_hc3) = 0.
1223         f0(p_hc5) = 0.
1224         f0(p_hc8) = 0.
1225         f0(p_olt) = 0.
1226         f0(p_oli) = 0.
1227         f0(p_tol) = 0.
1228         f0(p_csl) = 0.
1229         f0(p_xyl) = 0.
1230         f0(p_iso) = 0.
1231         f0(p_hno3) = 0.
1232         f0(p_ora1) = 0.
1233         f0(p_ora2) = 0.
1234         f0(p_nh3) = 0.
1235         f0(p_n2o5) = 1.
1236         if(p_ol2.gt.1)f0(p_ol2) = 0.
1237         f0(p_par) = 0.   !wig, 1-May-2007: for CBMZ
1238 !     DIFFUSION COEFFICIENTS
1239 !     [DV]=cm2/s (assumed: 1/SQRT(molar mass) when not known)
1240 
1241         dvj(p_no2) = 0.147
1242         dvj(p_no) = 0.183
1243         dvj(p_pan) = 0.091
1244         dvj(p_o3) = 0.175
1245         dvj(p_hcho) = 0.183
1246         dvj(p_aco3) = 0.115
1247         dvj(p_tpan) = 0.082
1248         dvj(p_hono) = 0.153
1249         dvj(p_no3) = 0.127
1250         dvj(p_hno4) = 0.113
1251         dvj(p_h2o2) = 0.171
1252         dvj(p_co) = 0.189
1253         dvj(p_ald) = 0.151
1254         dvj(p_op1) = 0.144
1255         dvj(p_op2) = 0.127
1256         dvj(p_paa) = 0.115
1257         dvj(p_ket) = 0.118
1258         dvj(p_gly) = 0.131
1259         dvj(p_mgly) = 0.118
1260         dvj(p_dcb) = 0.107
1261         dvj(p_onit) = 0.092
1262         dvj(p_so2) = 0.126
1263         dvj(p_eth) = 0.183
1264         dvj(p_hc3) = 0.151
1265         dvj(p_hc5) = 0.118
1266         dvj(p_hc8) = 0.094
1267         dvj(p_olt) = 0.154
1268         dvj(p_oli) = 0.121
1269         dvj(p_tol) = 0.104
1270         dvj(p_csl) = 0.096
1271         dvj(p_xyl) = 0.097
1272         dvj(p_iso) = 0.121
1273         dvj(p_hno3) = 0.126
1274         dvj(p_ora1) = 0.153
1275         dvj(p_ora2) = 0.124
1276         dvj(p_nh3) = 0.227
1277         dvj(p_n2o5) = 0.110
1278         dvj(p_ho) = 0.243
1279         dvj(p_ho2) = 0.174
1280         if(p_ol2.gt.1)dvj(p_ol2) = 0.189
1281         dvj(p_par) = 0.118   !wig, 1-May-2007: for CBMZ
1282         DO l = 1, numgas
1283           hstar4(l) = hstar(l) ! preliminary              
1284 ! Correction of diff. coeff
1285           dvj(l) = dvj(l)*(293.15/298.15)**1.75
1286           sc = 0.15/dvj(l) ! Schmidt Number at 20°C   
1287           dratio(l) = 0.242/dvj(l) !                                            ! of water vapor and gas at
1288 ! Ratio of diffusion coeffi
1289           scpr23(l) = (sc/0.72)**(2./3.) ! (Schmidt # / Prandtl #)**
1290         END DO
1291 
1292 
1293 
1294 !     DATA FOR AEROSOL PARTICLE DEPOSITION FOR THE MODEL OF
1295 !     J. W. ERISMAN, A. VAN PUL AND P. WYERS
1296 !     ATMOSPHERIC ENVIRONMENT 28 (1994), 2595-2607
1297 
1298 !     vd = (u* / k) * CORRECTION FACTORS
1299 
1300 !     CONSTANT K FOR LANDUSE TYPES:
1301 ! urban and built-up land                  
1302         kpart(1) = 500.
1303 ! dryland cropland and pasture             
1304         kpart(2) = 500.
1305 ! irrigated cropland and pasture           
1306         kpart(3) = 500.
1307 ! mixed dryland/irrigated cropland and past
1308         kpart(4) = 500.
1309 ! cropland/grassland mosaic                
1310         kpart(5) = 500.
1311 ! cropland/woodland mosaic                 
1312         kpart(6) = 100.
1313 ! grassland                                
1314         kpart(7) = 500.
1315 ! shrubland                                
1316         kpart(8) = 500.
1317 ! mixed shrubland/grassland                
1318         kpart(9) = 500.
1319 ! savanna                                  
1320         kpart(10) = 500.
1321 ! deciduous broadleaf forest               
1322         kpart(11) = 100.
1323 ! deciduous needleleaf forest              
1324         kpart(12) = 100.
1325 ! evergreen broadleaf forest               
1326         kpart(13) = 100.
1327 ! evergreen needleleaf forest              
1328         kpart(14) = 100.
1329 ! mixed forest                             
1330         kpart(15) = 100.
1331 ! water bodies                             
1332         kpart(16) = 500.
1333 ! herbaceous wetland                       
1334         kpart(17) = 500.
1335 ! wooded wetland                           
1336         kpart(18) = 500.
1337 ! barren or sparsely vegetated             
1338         kpart(19) = 500.
1339 ! herbaceous tundra                        
1340         kpart(20) = 500.
1341 ! wooded tundra                            
1342         kpart(21) = 100.
1343 ! mixed tundra                             
1344         kpart(22) = 500.
1345 ! bare ground tundra                       
1346         kpart(23) = 500.
1347 ! snow or ice                              
1348         kpart(24) = 500.
1349 !     Comments:
1350         kpart(25) = 500.
1351 !     Erisman et al. (1994) give
1352 !     k = 500 for low vegetation and k = 100 for forests.
1353 
1354 !     For desert k = 500 is taken according to measurements
1355 !     on bare soil by
1356 !     J. Fontan, A. Lopez, E. Lamaud and A. Druilhet (1997)
1357 !     Vertical Flux Measurements of the Submicronic Aerosol Particles
1358 !     and Parametrisation of the Dry Deposition Velocity
1359 !     in: Biosphere-Atmosphere Exchange of Pollutants
1360 !     and Trace Substances
1361 !     Editor: S. Slanina. Springer-Verlag Berlin, Heidelberg, 1997
1362 !     pp. 381-390
1363 
1364 !     For coniferous forest the Erisman value of  k = 100 is taken.
1365 !     Measurements of Erisman et al. (1997) in a coniferous forest
1366 !     in the Netherlands, lead to values of k between 20 and 38
1367 !     (Atmospheric Environment 31 (1997), 321-332).
1368 !     However, these high values of vd may be reached during
1369 !     instable cases. The eddy correlation measurements
1370 !     of Gallagher et al. (1997) made during the same experiment
1371 !     show for stable cases (L>0) values of k between 200 and 250
1372 !     at minimum (Atmospheric Environment 31 (1997), 359-373).
1373 !     Fontan et al. (1997) found k = 250 in a forest
1374 !     of maritime pine in southwestern France.
1375 
1376 !     For gras, model calculations of Davidson et al. support
1377 !     the value of 500.
1378 !     C. I. Davidson, J. M. Miller and M. A. Pleskov
1379 !     The Influence of Surface Structure on Predicted Particles
1380 !     Dry Deposition to Natural Gras Canopies
1381 !     Water, Air, and Soil Pollution 18 (1982) 25-43
1382 
1383 !     Snow covered surface: The experiment of Ibrahim et al. (1983)
1384 !     gives k = 436 for 0.7 um diameter particles.
1385 !     The deposition velocity of Milford and Davidson (1987)
1386 !     gives k = 154 for continental sulfate aerosol.
1387 !     M. Ibrahim, L. A. Barrie and F. Fanaki
1388 !     Atmospheric Environment 17 (1983), 781-788
1389 
1390 !     J. B. Milford and C. I. Davidson
1391 !     The Sizes of Particulate Sulfate and Nitrate in the Atmosphere
1392 !     - A Review
1393 !     JAPCA 37 (1987), 125-134
1394 ! no data                                  
1395 !       WRITE (0,*) ' return from rcread '
1396 !     *********************************************************
1397 
1398 !     Simplified landuse scheme for deposition and biogenic emission
1399 !     subroutines
1400 !     (ISWATER and ISICE are already defined elsewhere,
1401 !     therefore water and ice are not considered here)
1402 
1403 !     1 urban or bare soil
1404 !     2 agricultural
1405 !     3 grassland
1406 !     4 deciduous forest
1407 !     5 coniferous and mixed forest
1408 !     6 other natural landuse categories
1409 
1410 
1411         IF (mminlu=='OLD ') THEN
1412           ixxxlu(1) = 1
1413           ixxxlu(2) = 2
1414           ixxxlu(3) = 3
1415           ixxxlu(4) = 4
1416           ixxxlu(5) = 5
1417           ixxxlu(6) = 5
1418           ixxxlu(7) = 0
1419           ixxxlu(8) = 6
1420           ixxxlu(9) = 1
1421           ixxxlu(10) = 6
1422           ixxxlu(11) = 0
1423           ixxxlu(12) = 4
1424           ixxxlu(13) = 6
1425         END IF
1426         IF (mminlu=='USGS') THEN
1427           ixxxlu(1) = 1
1428           ixxxlu(2) = 2
1429           ixxxlu(3) = 2
1430           ixxxlu(4) = 2
1431           ixxxlu(5) = 2
1432           ixxxlu(6) = 4
1433           ixxxlu(7) = 3
1434           ixxxlu(8) = 6
1435           ixxxlu(9) = 3
1436           ixxxlu(10) = 6
1437           ixxxlu(11) = 4
1438           ixxxlu(12) = 5
1439           ixxxlu(13) = 4
1440           ixxxlu(14) = 5
1441           ixxxlu(15) = 5
1442           ixxxlu(16) = 0
1443           ixxxlu(17) = 6
1444           ixxxlu(18) = 4
1445           ixxxlu(19) = 1
1446           ixxxlu(20) = 6
1447           ixxxlu(21) = 4
1448           ixxxlu(22) = 6
1449           ixxxlu(23) = 1
1450           ixxxlu(24) = 0
1451           ixxxlu(25) = 1
1452         END IF
1453         IF (mminlu=='SiB ') THEN
1454           ixxxlu(1) = 4
1455           ixxxlu(2) = 4
1456           ixxxlu(3) = 4
1457           ixxxlu(4) = 5
1458           ixxxlu(5) = 5
1459           ixxxlu(6) = 6
1460           ixxxlu(7) = 3
1461           ixxxlu(8) = 6
1462           ixxxlu(9) = 6
1463           ixxxlu(10) = 6
1464           ixxxlu(11) = 1
1465           ixxxlu(12) = 2
1466           ixxxlu(13) = 6
1467           ixxxlu(14) = 1
1468           ixxxlu(15) = 0
1469           ixxxlu(16) = 0
1470           ixxxlu(17) = 1
1471         END IF
1472         RETURN
1473       END SUBROUTINE dep_init
1474     END MODULE module_dep_simple