Line data Source code
1 : !----------------------------------------------------------------------------------
2 : ! low level utility module for cloud aerosols
3 : !
4 : ! Created by Francis Vitt
5 : !----------------------------------------------------------------------------------
6 : module cldaero_mod
7 :
8 : use shr_kind_mod, only : r8 => shr_kind_r8
9 : use ppgrid, only : pcols, pver
10 :
11 : implicit none
12 : private
13 :
14 : public :: cldaero_uptakerate
15 : public :: cldaero_conc_t
16 : public :: cldaero_allocate
17 : public :: cldaero_deallocate
18 :
19 : type cldaero_conc_t
20 : real(r8), pointer :: so4c(:,:)
21 : real(r8), pointer :: nh4c(:,:)
22 : real(r8), pointer :: no3c(:,:)
23 : real(r8), pointer :: xlwc(:,:)
24 : real(r8) :: so4_fact
25 : end type cldaero_conc_t
26 :
27 : contains
28 :
29 : !----------------------------------------------------------------------------------
30 : !----------------------------------------------------------------------------------
31 0 : function cldaero_allocate( ) result( cldconc )
32 : type(cldaero_conc_t), pointer:: cldconc
33 :
34 0 : allocate( cldconc )
35 0 : allocate( cldconc%so4c(pcols,pver) )
36 0 : allocate( cldconc%nh4c(pcols,pver) )
37 0 : allocate( cldconc%no3c(pcols,pver) )
38 0 : allocate( cldconc%xlwc(pcols,pver) )
39 :
40 0 : cldconc%so4c(:,:) = 0._r8
41 0 : cldconc%nh4c(:,:) = 0._r8
42 0 : cldconc%no3c(:,:) = 0._r8
43 0 : cldconc%xlwc(:,:) = 0._r8
44 0 : cldconc%so4_fact = 2._r8
45 :
46 0 : end function cldaero_allocate
47 :
48 : !----------------------------------------------------------------------------------
49 : !----------------------------------------------------------------------------------
50 0 : subroutine cldaero_deallocate( cldconc )
51 : type(cldaero_conc_t), pointer :: cldconc
52 :
53 0 : if ( associated(cldconc%so4c) ) then
54 0 : deallocate(cldconc%so4c)
55 0 : nullify(cldconc%so4c)
56 : endif
57 :
58 0 : if ( associated(cldconc%nh4c) ) then
59 0 : deallocate(cldconc%nh4c)
60 0 : nullify(cldconc%nh4c)
61 : endif
62 :
63 0 : if ( associated(cldconc%no3c) ) then
64 0 : deallocate(cldconc%no3c)
65 0 : nullify(cldconc%no3c)
66 : endif
67 :
68 0 : if ( associated(cldconc%xlwc) ) then
69 0 : deallocate(cldconc%xlwc)
70 0 : nullify(cldconc%xlwc)
71 : endif
72 :
73 0 : deallocate( cldconc )
74 : nullify( cldconc )
75 :
76 0 : end subroutine cldaero_deallocate
77 :
78 : !----------------------------------------------------------------------------------
79 : ! utility function for cloud-borne aerosols
80 : !----------------------------------------------------------------------------------
81 :
82 0 : function cldaero_uptakerate( xl, cldnum, cfact, cldfrc, tfld, press ) result( uptkrate )
83 : use mo_constants, only : pi
84 :
85 : real(r8), intent(in) :: xl, cldnum, cfact, cldfrc, tfld, press
86 :
87 : real(r8) :: uptkrate
88 :
89 : real(r8) :: &
90 : rad_cd, radxnum_cd, num_cd, &
91 : gasdiffus, gasspeed, knudsen, &
92 : fuchs_sutugin, volx34pi_cd
93 :
94 : !-----------------------------------------------------------------------
95 : ! compute uptake of h2so4 and msa to cloud water
96 : !
97 : ! first-order uptake rate is
98 : ! 4*pi*(drop radius)*(drop number conc)
99 : ! *(gas diffusivity)*(fuchs sutugin correction)
100 :
101 : ! num_cd = (drop number conc in 1/cm^3)
102 0 : num_cd = 1.0e-3_r8*cldnum*cfact/cldfrc
103 0 : num_cd = max( num_cd, 0.0_r8 )
104 :
105 : ! rad_cd = (drop radius in cm), computed from liquid water and drop number,
106 : ! then bounded by 0.5 and 50.0 micrometers
107 : ! radxnum_cd = (drop radius)*(drop number conc)
108 : ! volx34pi_cd = (3/4*pi) * (liquid water volume in cm^3/cm^3)
109 :
110 0 : volx34pi_cd = xl*0.75_r8/pi
111 :
112 : ! following holds because volx34pi_cd = num_cd*(rad_cd**3)
113 0 : radxnum_cd = (volx34pi_cd*num_cd*num_cd)**0.3333333_r8
114 :
115 : ! apply bounds to rad_cd to avoid the occasional unphysical value
116 0 : if (radxnum_cd .le. volx34pi_cd*4.0e4_r8) then
117 : radxnum_cd = volx34pi_cd*4.0e4_r8
118 : rad_cd = 50.0e-4_r8
119 0 : else if (radxnum_cd .ge. volx34pi_cd*4.0e8_r8) then
120 : radxnum_cd = volx34pi_cd*4.0e8_r8
121 : rad_cd = 0.5e-4_r8
122 : else
123 0 : rad_cd = radxnum_cd/num_cd
124 : end if
125 :
126 : ! gasdiffus = h2so4 gas diffusivity from mosaic code (cm^2/s)
127 : ! (pmid must be Pa)
128 0 : gasdiffus = 0.557_r8 * (tfld**1.75_r8) / press
129 :
130 : ! gasspeed = h2so4 gas mean molecular speed from mosaic code (cm/s)
131 0 : gasspeed = 1.455e4_r8 * sqrt(tfld/98.0_r8)
132 :
133 : ! knudsen number
134 0 : knudsen = 3.0_r8*gasdiffus/(gasspeed*rad_cd)
135 :
136 : ! following assumes accomodation coefficient = 0.65
137 : ! (Adams & Seinfeld, 2002, JGR, and references therein)
138 : ! fuchs_sutugin = (0.75*accom*(1. + knudsen)) /
139 : ! (knudsen*(1.0 + knudsen + 0.283*accom) + 0.75*accom)
140 : fuchs_sutugin = (0.4875_r8*(1._r8 + knudsen)) / &
141 0 : (knudsen*(1.184_r8 + knudsen) + 0.4875_r8)
142 :
143 : ! instantaneous uptake rate
144 0 : uptkrate = 12.56637_r8*radxnum_cd*gasdiffus*fuchs_sutugin
145 :
146 0 : end function cldaero_uptakerate
147 :
148 0 : end module cldaero_mod
|