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