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