1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
3 | |
---|
4 | subroutine get_settling(itime,xt,yt,zt,nsp,settling) |
---|
5 | ! i i i i i o |
---|
6 | !***************************************************************************** |
---|
7 | ! * |
---|
8 | ! This subroutine calculates particle settling velocity. * |
---|
9 | ! * |
---|
10 | ! Author: A. Stohl * |
---|
11 | ! * |
---|
12 | ! May 2010 * |
---|
13 | ! * |
---|
14 | ! Improvement over traditional settling calculation in FLEXPART: * |
---|
15 | ! generalize to higher Reynolds numbers and also take into account the * |
---|
16 | ! temperature dependence of dynamic viscosity. * |
---|
17 | ! * |
---|
18 | ! Based on: * |
---|
19 | ! Naeslund E., and Thaning, L. (1991): On the settling velocity in a * |
---|
20 | ! nonstationary atmosphere, Aerosol Science and Technology 14, 247-256. * |
---|
21 | ! * |
---|
22 | !***************************************************************************** |
---|
23 | ! * |
---|
24 | ! Variables: * |
---|
25 | ! itime [s] current temporal position * |
---|
26 | ! xt,yt,zt coordinates position for which wind data shall be cal- * |
---|
27 | ! culated * |
---|
28 | ! * |
---|
29 | ! Constants: * |
---|
30 | ! * |
---|
31 | !***************************************************************************** |
---|
32 | |
---|
33 | use par_mod |
---|
34 | use com_mod |
---|
35 | |
---|
36 | implicit none |
---|
37 | |
---|
38 | integer :: itime,indz |
---|
39 | real :: xt,yt,zt |
---|
40 | |
---|
41 | ! Auxiliary variables needed for interpolation |
---|
42 | real :: dz1,dz2,dz |
---|
43 | real :: rho1(2),tt1(2),temperature,airdens,vis_dyn,vis_kin,viscosity |
---|
44 | real :: settling,settling_old,reynolds,c_d |
---|
45 | integer :: i,n,nix,njy,indzh,nsp |
---|
46 | |
---|
47 | |
---|
48 | !***************************************************************************** |
---|
49 | ! 1. Interpolate temperature and density: nearest neighbor interpolation sufficient |
---|
50 | !***************************************************************************** |
---|
51 | |
---|
52 | nix=int(xt) |
---|
53 | njy=int(yt) |
---|
54 | |
---|
55 | ! Determine the level below the current position for u,v |
---|
56 | !******************************************************* |
---|
57 | |
---|
58 | do i=2,nz |
---|
59 | if (height(i).gt.zt) then |
---|
60 | indz=i-1 |
---|
61 | goto 6 |
---|
62 | endif |
---|
63 | end do |
---|
64 | 6 continue |
---|
65 | |
---|
66 | |
---|
67 | ! Vertical distance to the level below and above current position |
---|
68 | !**************************************************************** |
---|
69 | |
---|
70 | dz=1./(height(indz+1)-height(indz)) |
---|
71 | dz1=(zt-height(indz))*dz |
---|
72 | dz2=(height(indz+1)-zt)*dz |
---|
73 | |
---|
74 | |
---|
75 | ! Bilinear horizontal interpolation |
---|
76 | !********************************** |
---|
77 | |
---|
78 | ! Loop over 2 levels |
---|
79 | !******************* |
---|
80 | |
---|
81 | do n=1,2 |
---|
82 | indzh=indz+n-1 |
---|
83 | rho1(n)=rho(nix,njy,indzh,1) |
---|
84 | tt1(n)=tt(nix,njy,indzh,1) |
---|
85 | end do |
---|
86 | |
---|
87 | |
---|
88 | ! Linear vertical interpolation |
---|
89 | !****************************** |
---|
90 | |
---|
91 | temperature=dz2*tt1(1)+dz1*tt1(2) |
---|
92 | airdens=dz2*rho1(1)+dz1*rho1(2) |
---|
93 | |
---|
94 | |
---|
95 | vis_dyn=viscosity(temperature) |
---|
96 | vis_kin=vis_dyn/airdens |
---|
97 | |
---|
98 | reynolds=dquer(nsp)/1.e6*abs(vsetaver(nsp))/vis_kin |
---|
99 | |
---|
100 | |
---|
101 | |
---|
102 | ! Iteration to determine both Reynolds number and settling velocity |
---|
103 | !****************************************************************** |
---|
104 | |
---|
105 | settling_old=vsetaver(nsp) ! initialize iteration with Stokes' law, constant viscosity estimate |
---|
106 | |
---|
107 | do i=1,20 ! do a few iterations |
---|
108 | |
---|
109 | if (reynolds.lt.1.917) then |
---|
110 | c_d=24./reynolds |
---|
111 | else if (reynolds.lt.500.) then |
---|
112 | c_d=18.5/(reynolds**0.6) |
---|
113 | else |
---|
114 | c_d=0.44 |
---|
115 | endif |
---|
116 | |
---|
117 | settling=-1.* & |
---|
118 | sqrt(4*ga*dquer(nsp)/1.e6*density(nsp)*cunningham(nsp)/ & |
---|
119 | (3.*c_d*airdens)) |
---|
120 | |
---|
121 | if (abs((settling-settling_old)/settling).lt.0.01) goto 11 ! stop iteration |
---|
122 | |
---|
123 | reynolds=dquer(nsp)/1.e6*abs(settling)/vis_kin |
---|
124 | settling_old=settling |
---|
125 | end do |
---|
126 | |
---|
127 | 11 continue |
---|
128 | |
---|
129 | end subroutine get_settling |
---|