Line data Source code
1 : ! Include shortname defintions, so that the F77 code does not have to be modified to
2 : ! reference the CARMA structure.
3 : #include "carma_globaer.h"
4 :
5 :
6 : !!-----------------------------------------------------------------------
7 : !!
8 : !! Purpose: Calculating the surface resistance, using the PBL parameter.
9 : !!
10 : !! Method: Zhang(2001), Atmospheric Environment
11 : !!
12 : !!
13 : !! @author Tianyi Fan
14 : !! @version Nov-2010
15 : !!
16 :
17 175239400 : subroutine calcrs(carma, cstate, ustar, tmp, radi, cc, vfall, rs, landidx, rc)
18 : use carma_precision_mod
19 : use carmastate_mod
20 : use carma_enums_mod
21 : use carma_types_mod
22 : use carma_mod
23 : use carma_constants_mod, only: BK, PI, GRAV
24 : !-----------------------------------------------------------------------
25 : implicit none
26 : !-----------------------------------------------------------------------
27 : !
28 : ! Arguments:
29 : !
30 : type(carma_type), intent(in) :: carma !! the carma object
31 : type(carmastate_type), intent(in) :: cstate !! the carma state object
32 : real(kind=f), intent(in) :: ustar !! friction velocity [cm/s]
33 : real(kind=f), intent(in) :: tmp !! temperature [K]
34 : real(kind=f), intent(in) :: radi !! radius of the constitutent [cm]
35 : real(kind=f), intent(in) :: cc !! slip correction factor
36 : real(kind=f), intent(in) :: vfall !! gravitational settling velocity, [cm/s]
37 : real(kind=f), intent(out) :: rs !! surface resistance [s/cm]
38 : integer, intent(in) :: landidx !! landscape index, 1=land, 2=ocean, 3=sea ice
39 : integer, intent(inout) :: rc !! return code, negative indicates failure
40 :
41 :
42 :
43 : ! Local variables
44 : real(kind=f) :: ebrn ! Brownian diffusion collection efficiency
45 : real(kind=f) :: eimp ! Impaction collection efficiency
46 : real(kind=f) :: eint ! Interception collection efficiency
47 : real(kind=f) :: db ! Brownian diffusivity
48 : real(kind=f) :: sc ! Schmidt number
49 : real(kind=f) :: st ! Stokes number
50 : real(kind=f) :: rhoadry ! dry air density [g/cm3]
51 : real(kind=f) :: eta ! kinematic viscosity of air [cm2/s]
52 : real(kind=f), parameter :: xkar = 0.4_f ! Von Karman's constant
53 : real(kind=f), parameter :: eps0 = 3._f ! empirical constant for rs, 3.0 in [Zhang, 2001], 1.0 in [Seinfeld and Pandis]
54 :
55 : ! exponent in the eb dependence on sc, 2/3 in [Seinfeld and Pandis, 1998], 1/2 in [Lewis and Schwartz, 2004]
56 : real(kind=f) :: lam
57 :
58 : integer :: ibot
59 :
60 175239400 : if (igridv .eq. I_CART) then
61 : ibot = 1
62 : else
63 175239400 : ibot = NZ
64 : end if
65 :
66 : ! Unit conversion
67 175239400 : rhoadry = rhoa(ibot) / zmet(ibot) ! [g/cm3]
68 175239400 : eta = rmu(ibot) / rhoadry ! rmu, aerodynamic viscosity of air [g/cm/s]
69 :
70 175239400 : if (landidx .eq. 1) then
71 : lam = 2._f / 3._f
72 : else
73 119344820 : lam = 1._f / 2._f
74 : end if
75 :
76 : ! Surface Resistance = Brownian + Impaction + Interception
77 :
78 : ! ** Brownian diffusion
79 175239400 : db = BK * tmp * cc / (6._f * PI * rmu(ibot) * radi) ! [cm2/s]
80 :
81 175239400 : sc = eta / db ! [-]
82 175239400 : ebrn = sc**(-lam)
83 :
84 : ! ** Impaction
85 175239400 : st = vfall * ustar**2 / (GRAV * eta) ! [-]
86 :
87 : ! [Slinn, 1982]
88 : ! eimp = 10. ** (-3._f/st)
89 :
90 : ! [Peters and Eiden, 1992]
91 175239400 : eimp = (st / (0.8_f + st))**2
92 : ! eimp = max(eimp, 1.e-10_f)
93 :
94 : ! ** Interception
95 : !
96 : ! NOTE: Interception is not currently considered for ocean and ice.
97 175239400 : if (landidx .eq. 1) then
98 : ! eint = 0.3_f * (0.01_f * radi * 1.e-2_f / (radi * 1.e-2_f + 1.e-5_f) + 0.99_f *radi*1.e-2_f / (radi*1.e-2_f + 8.e-4_f))
99 55894580 : eint = 0.3_f * (0.01_f * radi / (radi + 1.e-3_f) + 0.99_f * radi / (radi + 8.e-2_f))
100 : else
101 : eint = 0._f
102 : end if
103 :
104 175239400 : if (ustar > 0._f) then
105 175239400 : rs = 1._f / (eps0 * ustar * (ebrn + eimp + eint )) ! [s/cm]
106 : else
107 0 : rs = 0._f
108 : end if
109 :
110 175239400 : return
111 175239400 : end subroutine calcrs
|