source: flexpart.git/src/partdep.f90

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

add SPDX-License-Identifier to all .f90 files

  • Property mode set to 100644
File size: 5.4 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine partdep(nc,density,fract,schmi,vset,ra,ustar,nyl,vdep)
5  !                   i     i      i     i    i   i    i    i  i/o
6  !*****************************************************************************
7  !                                                                            *
8  !      Calculation of the dry deposition velocities of particles.            *
9  !      This routine is based on Stokes' law for considering settling and     *
10  !      assumes constant dynamic viscosity of the air.                        *
11  !                                                                            *
12  !     AUTHOR: Andreas Stohl, 12 November 1993                                *
13  !                            Update: 20 December 1996                        *
14  !                                                                            *
15  !     Literature:                                                            *
16  !     [1]  Hicks/Baldocchi/Meyers/Hosker/Matt (1987), A Preliminary          *
17  !             Multiple Resistance Routine for Deriving Dry Deposition        *
18  !             Velocities from Measured Quantities.                           *
19  !             Water, Air and Soil Pollution 36 (1987), pp.311-330.           *
20  !     [2]  Slinn (1982), Predictions for Particle Deposition to              *
21  !             Vegetative Canopies. Atm.Env.16-7 (1982), pp.1785-1794.        *
22  !     [3]  Slinn/Slinn (1980),  Predictions for Particle Deposition on       *
23  !             Natural Waters. Atm.Env.14 (1980), pp.1013-1016.               *
24  !     [4]  Scire/Yamartino/Carmichael/Chang (1989),                          *
25  !             CALGRID: A Mesoscale Photochemical Grid Model.                 *
26  !             Vol II: User's Guide. (Report No.A049-1, June, 1989)           *
27  !     [5]  Langer M. (1992): Ein einfaches Modell zur Abschaetzung der       *
28  !             Depositionsgeschwindigkeit von Teilchen und Gasen.             *
29  !             Internal report.                                               *
30  !                                                                            *
31  !*****************************************************************************
32  !                                                                            *
33  ! Variables:                                                                 *
34  ! alpha                help variable                                         *
35  ! fract(nc,ni)         mass fraction of each diameter interval               *
36  ! lpdep(nc)            1 for particle deposition, 0 else                     *
37  ! nc                   actual number of chemical components                  *
38  ! ni                   number of diameter intervals, for which vdepj is calc.*
39  ! rdp [s/m]            deposition layer resistance                           *
40  ! ra [s/m]             aerodynamical resistance                              *
41  ! schmi(nc,ni)         Schmidt number**2/3 of each diameter interval         *
42  ! stokes               Stokes number                                         *
43  ! ustar [m/s]          friction velocity                                     *
44  ! vdep(nc) [m/s]       deposition velocities of all components               *
45  ! vdepj [m/s]          help, deposition velocity of 1 interval               *
46  ! vset(nc,ni)          gravitational settling velocity of each interval      *
47  !                                                                            *
48  ! Constants:                                                                 *
49  ! nc                   number of chemical species                            *
50  ! ni                   number of diameter intervals, for which deposition    *
51  !                      is calculated                                         *
52  !                                                                            *
53  !*****************************************************************************
54
55  use par_mod
56  use com_mod, only: debug_mode
57
58  implicit none
59
60  real :: density(maxspec),schmi(maxspec,ni),fract(maxspec,ni)
61  real :: vset(maxspec,ni)
62  real :: vdep(maxspec),stokes,vdepj,rdp,ustar,alpha,ra,nyl
63  real,parameter :: eps=1.e-5
64  integer :: ic,j,nc
65
66
67  do ic=1,nc                  ! loop over all species
68    if (density(ic).gt.0.) then
69      do j=1,ni              ! loop over all diameter intervals
70        if (ustar.gt.eps) then
71
72  ! Stokes number for each diameter interval
73  !*****************************************
74
75          stokes=vset(ic,j)/ga*ustar*ustar/nyl
76          alpha=-3./stokes
77
78  ! Deposition layer resistance
79  !****************************
80
81          if (alpha.le.log10(eps)) then
82            rdp=1./(schmi(ic,j)*ustar)
83          else
84            rdp=1./((schmi(ic,j)+10.**alpha)*ustar)
85          endif
86          vdepj=vset(ic,j)+1./(ra+rdp+ra*rdp*vset(ic,j))
87        else
88          vdepj=vset(ic,j)
89        endif
90
91  ! deposition velocities of each interval are weighted with mass fraction
92  !***********************************************************************
93
94        vdep(ic)=vdep(ic)+vdepj*fract(ic,j)
95       
96        !print*, 'partdep:113: ic', ic, 'vdep', vdep
97        !stop
98      end do
99    endif
100
101
102  end do
103
104  !if (debug_mode) then
105  !  print*, 'partdep:122:'
106  !  write(*,*) (vdep(ic), ic=1,nc)
107    !stop
108  !endif
109
110end subroutine partdep
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG