module_dep_simple.F

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