source: flexpart.git/src/getvdep_nests.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 5.9 KB
Line 
1subroutine getvdep_nests(n,ix,jy,ust,temp,pa, &
2       L,gr,rh,rr,snow,vdepo,lnest)
3  !                   i i  i   i   i   i  i i  i  i    i o i
4  !*****************************************************************************
5  !                                                                            *
6  !  This routine calculates the dry deposition velocities.                    *
7  !                                                                            *
8  !     Author: A. Stohl                                                       *
9  !                                                                            *
10  !     20 December 1996                                                       *
11  !     Sabine Eckhardt, Jan 07                                                *
12  !     if the latitude is negative: add half a year to the julian day         *
13  !                                                                            *
14  !*****************************************************************************
15  !                                                                            *
16  ! Variables:                                                                 *
17  ! gr [W/m2]         global radiation                                         *
18  ! L [m]             Obukhov length                                           *
19  ! nyl               kinematic viscosity                                      *
20  ! pa [Pa]           surface air pressure                                     *
21  ! ra [s/m]          aerodynamic resistance                                   *
22  ! raquer [s/m]      average aerodynamic resistance                           *
23  ! rh [0-1]          relative humidity                                        *
24  ! rhoa              density of the air                                       *
25  ! rr [mm/h]         precipitation rate                                       *
26  ! temp [K]          2m temperature                                           *
27  ! tc [C]            2m temperature                                           *
28  ! ust [m/s]         friction velocity                                        *
29  ! snow [m of water equivalent] snow depth                                    *
30  ! xlanduse          fractions of numclasS landuses for each model grid point *
31  !                                                                            *
32  !*****************************************************************************
33
34  use par_mod
35  use com_mod
36
37  implicit none
38
39  integer :: yyyymmdd,hhmmss,yyyy,mmdd,n,lseason,i,j,ix,jy,lnest
40  real :: vdepo(maxspec),vd,rb(maxspec),rc(maxspec),raquer,ylat
41  real :: raerod,ra,ust,temp,tc,pa,L,gr,rh,rr,myl,nyl,rhoa,diffh2o,snow
42  real :: slanduse(numclass)
43  real,parameter :: eps=1.e-5
44  real(kind=dp) :: jul
45
46  ! Calculate month and determine the seasonal category
47  !****************************************************
48
49  jul=bdate+real(wftime(n),kind=dp)/86400._dp
50
51  ylat=jy*dy+ylat0
52  if (ylat.lt.0) then
53      jul=jul+365/2
54  endif
55
56
57  call caldate(jul,yyyymmdd,hhmmss)
58  yyyy=yyyymmdd/10000
59  mmdd=yyyymmdd-10000*yyyy
60
61  if ((ylat.gt.-20).and.(ylat.lt.20)) then
62     mmdd=600 ! summer
63  endif
64
65  if ((mmdd.ge.1201).or.(mmdd.le.301)) then
66    lseason=4
67  else if ((mmdd.ge.1101).or.(mmdd.le.331)) then
68    lseason=3
69  else if ((mmdd.ge.401).and.(mmdd.le.515)) then
70    lseason=5
71  else if ((mmdd.ge.516).and.(mmdd.le.915)) then
72    lseason=1
73  else
74    lseason=2
75  endif
76
77  ! Calculate diffusivity of water vapor
78  !************************************
79  diffh2o=2.11e-5*(temp/273.15)**1.94*(101325/pa)
80
81  ! Conversion of temperature from K to C
82  !**************************************
83
84  tc=temp-273.15
85
86  ! Calculate dynamic viscosity
87  !****************************
88
89  if (tc.lt.0) then
90    myl=(1.718+0.0049*tc-1.2e-05*tc**2)*1.e-05
91  else
92    myl=(1.718+0.0049*tc)*1.e-05
93  endif
94
95  ! Calculate kinematic viscosity
96  !******************************
97
98  rhoa=pa/(287.*temp)
99  nyl=myl/rhoa
100
101
102  ! 0. Set all deposition velocities zero
103  !**************************************
104
105  do i=1,nspec
106    vdepo(i)=0.
107  end do
108
109
110  ! 1. Compute surface layer resistances rb
111  !****************************************
112
113  call getrb(nspec,ust,nyl,diffh2o,reldiff,rb)
114
115  ! change for snow
116  do j=1,numclass
117    if (snow.gt.0.001) then ! 10 mm
118       if (j.eq.12) then
119         slanduse(j)=1.
120       else
121         slanduse(j)=0.
122       endif
123    else
124       slanduse(j)=xlandusen(ix,jy,j,lnest)
125    endif
126  end do
127
128  raquer=0.
129  do j=1,numclass            ! loop over all landuse classes
130
131    if (slanduse(j).gt.eps)  then
132
133  ! 2. Calculate aerodynamic resistance ra
134  !***************************************
135
136      ra=raerod(L,ust,z0(j))
137      raquer=raquer+ra*slanduse(j)
138
139  ! 3. Calculate surface resistance for gases
140  !******************************************
141
142      call getrc(nspec,lseason,j,tc,gr,rh,rr,rc)
143
144  ! 4. Calculate deposition velocities for gases and ...
145  ! 5. ... sum deposition velocities for all landuse classes
146  !*********************************************************
147
148      do i=1,nspec
149        if (reldiff(i).gt.0.) then
150          if ((ra+rb(i)+rc(i)).gt.0.) then
151            vd=1./(ra+rb(i)+rc(i))
152  ! XXXXXXXXXXXXXXXXXXXXXXXXXX TEST
153  !         vd=1./rc(i)
154  ! XXXXXXXXXXXXXXXXXXXXXXXXXX TEST
155          else
156            vd=9.999
157          endif
158          vdepo(i)=vdepo(i)+vd*slanduse(j)
159        endif
160      end do
161    endif
162  end do
163
164
165  ! 6. Calculate deposition velocities for particles
166  !*************************************************
167
168  call partdep(nspec,density,fract,schmi,vset,raquer,ust,nyl,vdepo)
169
170
171  ! 7. If no detailed parameterization available, take constant deposition
172  !    velocity if that is available
173  !***********************************************************************
174
175  do i=1,nspec
176    if ((reldiff(i).lt.0.).and.(density(i).lt.0.).and. &
177         (dryvel(i).gt.0.)) then
178      vdepo(i)=dryvel(i)
179    endif
180  end do
181
182
183end subroutine getvdep_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG