Line data Source code
1 : module aero_deposition_cam
2 : !------------------------------------------------------------------------------
3 : ! Purpose:
4 : !
5 : ! Partition the contributions from aerosols of wet and dry
6 : ! deposition at the surface into the fields passed to the coupler.
7 : !------------------------------------------------------------------------------
8 :
9 : use shr_kind_mod, only: r8 => shr_kind_r8
10 : use shr_kind_mod, only: cl => shr_kind_cl
11 : use constituents, only: cnst_get_ind, pcnst
12 : use camsrfexch, only: cam_out_t
13 : use cam_abortutils,only: endrun
14 : use aerosol_properties_mod, only: aero_name_len
15 : use aerosol_properties_mod, only: aerosol_properties
16 :
17 : implicit none
18 :
19 : private
20 :
21 : ! Public interfaces
22 :
23 : public :: aero_deposition_cam_init
24 : public :: aero_deposition_cam_setwet
25 : public :: aero_deposition_cam_setdry
26 :
27 : ! Private module data
28 :
29 : integer :: bcphi_ndx( pcnst ) = -1
30 : integer :: bcphi_cnt = 0
31 : integer :: bcpho_ndx( pcnst ) = -1
32 : integer :: bcpho_cnt = 0
33 : integer :: ocphi_ndx( pcnst ) = -1
34 : integer :: ocphi_cnt = 0
35 : integer :: ocpho_ndx( pcnst ) = -1
36 : integer :: ocpho_cnt = 0
37 :
38 : class(aerosol_properties), pointer :: aero_props=>null()
39 : integer :: nele_tot=0 ! total number of aerosol elements
40 :
41 : ! bulk dust bins (meters)
42 :
43 : integer, parameter :: n_bulk_dst_bins = 4
44 :
45 : ! CAM4 bulk dust bin sizes (https://doi.org/10.1002/2013MS000279)
46 : real(r8), parameter :: bulk_dst_edges(n_bulk_dst_bins+1) = &
47 : (/0.1e-6_r8, 1.0e-6_r8, 2.5e-6_r8, 5.0e-6_r8, 10.e-6_r8/)
48 :
49 : contains
50 :
51 : !============================================================================
52 0 : subroutine aero_deposition_cam_init(aero_props_in)
53 :
54 : class(aerosol_properties),target, intent(in) :: aero_props_in
55 :
56 : integer :: pcnt, scnt
57 : character(len=*), parameter :: subrname = 'aero_deposition_cam_init'
58 :
59 : ! construct the aerosol properties object
60 0 : aero_props => aero_props_in
61 :
62 : ! set the cam constituent indices and determine the counts
63 : ! for the specified aerosol types
64 :
65 : ! black carbons
66 0 : call get_indices( type='black-c', hydrophilic=.true., indices=bcphi_ndx, count=bcphi_cnt )
67 0 : call get_indices( type='black-c', hydrophilic=.false., indices=bcpho_ndx, count=bcpho_cnt )
68 :
69 : ! primary and secondary organics
70 0 : call get_indices( type='p-organic',hydrophilic=.true., indices=ocphi_ndx, count=pcnt )
71 0 : call get_indices( type='s-organic',hydrophilic=.true., indices=ocphi_ndx(pcnt+1:), count=scnt )
72 0 : ocphi_cnt = pcnt+scnt
73 :
74 0 : call get_indices( type='p-organic',hydrophilic=.false., indices=ocpho_ndx, count=pcnt )
75 0 : call get_indices( type='s-organic',hydrophilic=.false., indices=ocpho_ndx(pcnt+1:), count=scnt )
76 0 : ocpho_cnt = pcnt+scnt
77 :
78 : ! total number of aerosol elements
79 0 : nele_tot = aero_props%ncnst_tot()
80 :
81 : contains
82 :
83 : !==========================================================================
84 : ! returns CAM constituent indices of the aerosol tracers (and count)
85 : !==========================================================================
86 0 : subroutine get_indices( type, hydrophilic, indices, count)
87 :
88 : character(len=*), intent(in) :: type
89 : logical, intent(in ) :: hydrophilic
90 : integer, intent(out) :: indices(:)
91 : integer, intent(out) :: count
92 :
93 : integer :: ibin,ispc, ndx, nspec
94 : character(len=aero_name_len) :: spec_type, spec_name
95 :
96 0 : count = 0
97 0 : indices(:) = -1
98 :
99 : ! loop through aerosol bins / modes
100 0 : do ibin = 1, aero_props%nbins()
101 :
102 : ! check if the bin/mode is hydrophilic
103 0 : if ( aero_props%hydrophilic(ibin) .eqv. hydrophilic ) then
104 0 : do ispc = 1, aero_props%nspecies(ibin)
105 :
106 0 : call aero_props%get(ibin,ispc, spectype=spec_type, specname=spec_name)
107 :
108 0 : if (spec_type==type) then
109 :
110 : ! get CAM constituent index
111 0 : call cnst_get_ind(spec_name, ndx, abort=.false.)
112 0 : if (ndx>0) then
113 0 : count = count+1
114 0 : indices(count) = ndx
115 : endif
116 :
117 : endif
118 :
119 : enddo
120 : endif
121 :
122 : enddo
123 :
124 0 : end subroutine get_indices
125 :
126 : end subroutine aero_deposition_cam_init
127 :
128 : !============================================================================
129 : ! Set surface wet deposition fluxes passed to coupler.
130 : !============================================================================
131 0 : subroutine aero_deposition_cam_setwet(aerdepwetis, aerdepwetcw, cam_out)
132 :
133 : ! Arguments:
134 : real(r8), intent(in) :: aerdepwetis(:,:) ! aerosol wet deposition (interstitial)
135 : real(r8), intent(in) :: aerdepwetcw(:,:) ! aerosol wet deposition (cloud water)
136 : type(cam_out_t), intent(inout) :: cam_out ! cam export state
137 :
138 : ! Local variables:
139 : integer :: i, ispec, ibin, mm, ndx
140 : integer :: ncol ! number of columns
141 :
142 0 : real(r8) :: dep_fluxes(nele_tot)
143 : real(r8) :: dst_fluxes(n_bulk_dst_bins)
144 : character(len=aero_name_len) :: specname, name_c
145 : integer :: errstat
146 : character(len=cl) :: errstr
147 :
148 0 : ncol = cam_out%ncol
149 :
150 0 : cam_out%bcphiwet(:) = 0._r8
151 0 : cam_out%ocphiwet(:) = 0._r8
152 0 : cam_out%dstwet1(:) = 0._r8
153 0 : cam_out%dstwet2(:) = 0._r8
154 0 : cam_out%dstwet3(:) = 0._r8
155 0 : cam_out%dstwet4(:) = 0._r8
156 :
157 : ! derive cam_out variables from deposition fluxes
158 : ! note: wet deposition fluxes are negative into surface,
159 : ! dry deposition fluxes are positive into surface.
160 : ! srf models want positive definite fluxes.
161 0 : do i = 1, ncol
162 :
163 : ! hydrophilic black carbon fluxes
164 0 : do ispec=1,bcphi_cnt
165 0 : cam_out%bcphiwet(i) = cam_out%bcphiwet(i) &
166 0 : - (aerdepwetis(i,bcphi_ndx(ispec))+aerdepwetcw(i,bcphi_ndx(ispec)))
167 : enddo
168 :
169 : ! hydrophobic black carbon fluxes
170 0 : do ispec=1,bcpho_cnt
171 0 : cam_out%bcphiwet(i) = cam_out%bcphiwet(i) &
172 0 : - (aerdepwetis(i,bcpho_ndx(ispec))+aerdepwetcw(i,bcpho_ndx(ispec)))
173 : enddo
174 :
175 : ! hydrophilic organic carbon fluxes
176 0 : do ispec=1,ocphi_cnt
177 0 : cam_out%ocphiwet(i) = cam_out%ocphiwet(i) &
178 0 : - (aerdepwetis(i,ocphi_ndx(ispec))+aerdepwetcw(i,ocphi_ndx(ispec)))
179 : enddo
180 :
181 : ! hydrophobic organic carbon fluxes
182 0 : do ispec=1,ocpho_cnt
183 0 : cam_out%ocphiwet(i) = cam_out%ocphiwet(i) &
184 0 : - (aerdepwetis(i,ocpho_ndx(ispec))+aerdepwetcw(i,ocpho_ndx(ispec)))
185 : enddo
186 :
187 : ! dust fluxes
188 :
189 0 : dep_fluxes = 0._r8
190 0 : dst_fluxes = 0._r8
191 :
192 0 : do ibin = 1,aero_props%nbins()
193 0 : do ispec = 0,aero_props%nmasses(ibin)
194 0 : if (ispec==0) then
195 0 : call aero_props%num_names(ibin, specname, name_c)
196 : else
197 0 : call aero_props%get(ibin,ispec, specname=specname)
198 : end if
199 0 : call cnst_get_ind(specname, ndx, abort=.false.)
200 0 : if (ndx>0) then
201 0 : mm = aero_props%indexer(ibin,ispec)
202 0 : dep_fluxes(mm) = - (aerdepwetis(i,ndx)+aerdepwetcw(i,ndx))
203 : end if
204 : end do
205 : end do
206 :
207 : ! rebin dust fluxes to bulk dust bins
208 0 : call aero_props%rebin_bulk_fluxes('dust', dep_fluxes, bulk_dst_edges, dst_fluxes, errstat, errstr)
209 0 : if (errstat/=0) then
210 0 : call endrun('aero_deposition_cam_setwet: '//trim(errstr))
211 : end if
212 :
213 0 : cam_out%dstwet1(i) = cam_out%dstwet1(i) + dst_fluxes(1)
214 0 : cam_out%dstwet2(i) = cam_out%dstwet2(i) + dst_fluxes(2)
215 0 : cam_out%dstwet3(i) = cam_out%dstwet3(i) + dst_fluxes(3)
216 0 : cam_out%dstwet4(i) = cam_out%dstwet4(i) + dst_fluxes(4)
217 :
218 : ! in rare cases, integrated deposition tendency is upward
219 0 : if (cam_out%bcphiwet(i) < 0._r8) cam_out%bcphiwet(i) = 0._r8
220 0 : if (cam_out%ocphiwet(i) < 0._r8) cam_out%ocphiwet(i) = 0._r8
221 0 : if (cam_out%dstwet1(i) < 0._r8) cam_out%dstwet1(i) = 0._r8
222 0 : if (cam_out%dstwet2(i) < 0._r8) cam_out%dstwet2(i) = 0._r8
223 0 : if (cam_out%dstwet3(i) < 0._r8) cam_out%dstwet3(i) = 0._r8
224 0 : if (cam_out%dstwet4(i) < 0._r8) cam_out%dstwet4(i) = 0._r8
225 :
226 : enddo
227 :
228 0 : end subroutine aero_deposition_cam_setwet
229 :
230 : !============================================================================
231 : ! Set surface dry deposition fluxes passed to coupler.
232 : !============================================================================
233 0 : subroutine aero_deposition_cam_setdry(aerdepdryis, aerdepdrycw, cam_out)
234 :
235 : ! Arguments:
236 : real(r8), intent(in) :: aerdepdryis(:,:) ! aerosol dry deposition (interstitial)
237 : real(r8), intent(in) :: aerdepdrycw(:,:) ! aerosol dry deposition (cloud water)
238 : type(cam_out_t), intent(inout) :: cam_out ! cam export state
239 :
240 : ! Local variables:
241 : integer :: i, ispec, ibin, mm, ndx
242 : integer :: ncol ! number of columns
243 :
244 0 : real(r8) :: dep_fluxes(nele_tot)
245 : real(r8) :: dst_fluxes(n_bulk_dst_bins)
246 : character(len=aero_name_len) :: specname, name_c
247 : integer :: errstat
248 : character(len=cl) :: errstr
249 :
250 0 : ncol = cam_out%ncol
251 :
252 0 : cam_out%bcphidry(:) = 0._r8
253 0 : cam_out%ocphidry(:) = 0._r8
254 0 : cam_out%bcphodry(:) = 0._r8
255 0 : cam_out%ocphodry(:) = 0._r8
256 0 : cam_out%dstdry1(:) = 0._r8
257 0 : cam_out%dstdry2(:) = 0._r8
258 0 : cam_out%dstdry3(:) = 0._r8
259 0 : cam_out%dstdry4(:) = 0._r8
260 :
261 : ! derive cam_out variables from deposition fluxes
262 : ! note: wet deposition fluxes are negative into surface,
263 : ! dry deposition fluxes are positive into surface.
264 : ! srf models want positive definite fluxes.
265 0 : do i = 1, ncol
266 :
267 : ! hydrophilic black carbon fluxes
268 0 : do ispec=1,bcphi_cnt
269 0 : cam_out%bcphidry(i) = cam_out%bcphidry(i) &
270 0 : + (aerdepdryis(i,bcphi_ndx(ispec))+aerdepdrycw(i,bcphi_ndx(ispec)))
271 : enddo
272 :
273 : ! hydrophobic black carbon fluxes
274 0 : do ispec=1,bcpho_cnt
275 0 : cam_out%bcphodry(i) = cam_out%bcphodry(i) &
276 0 : + (aerdepdryis(i,bcpho_ndx(ispec))+aerdepdrycw(i,bcpho_ndx(ispec)))
277 : enddo
278 :
279 : ! hydrophilic organic carbon fluxes
280 0 : do ispec=1,ocphi_cnt
281 0 : cam_out%ocphidry(i) = cam_out%ocphidry(i) &
282 0 : + (aerdepdryis(i,ocphi_ndx(ispec))+aerdepdrycw(i,ocphi_ndx(ispec)))
283 : enddo
284 :
285 : ! hydrophobic organic carbon fluxes
286 0 : do ispec=1,ocpho_cnt
287 0 : cam_out%ocphodry(i) = cam_out%ocphodry(i) &
288 0 : + (aerdepdryis(i,ocpho_ndx(ispec))+aerdepdrycw(i,ocpho_ndx(ispec)))
289 : enddo
290 :
291 : ! dust fluxes
292 :
293 0 : dep_fluxes = 0._r8
294 0 : dst_fluxes = 0._r8
295 :
296 0 : do ibin = 1,aero_props%nbins()
297 0 : do ispec = 0,aero_props%nspecies(ibin)
298 0 : if (ispec==0) then
299 0 : call aero_props%num_names(ibin, specname, name_c)
300 : else
301 0 : call aero_props%get(ibin,ispec, specname=specname)
302 : end if
303 0 : call cnst_get_ind(specname, ndx, abort=.false.)
304 0 : if (ndx>0) then
305 0 : mm = aero_props%indexer(ibin,ispec)
306 0 : dep_fluxes(mm) = aerdepdryis(i,ndx)+aerdepdrycw(i,ndx)
307 : end if
308 : end do
309 : end do
310 :
311 : ! rebin dust fluxes to bulk dust bins
312 0 : call aero_props%rebin_bulk_fluxes('dust', dep_fluxes, bulk_dst_edges, dst_fluxes, errstat, errstr)
313 0 : if (errstat/=0) then
314 0 : call endrun('aero_deposition_cam_setdry: '//trim(errstr))
315 : end if
316 :
317 0 : cam_out%dstdry1(i) = cam_out%dstdry1(i) + dst_fluxes(1)
318 0 : cam_out%dstdry2(i) = cam_out%dstdry2(i) + dst_fluxes(2)
319 0 : cam_out%dstdry3(i) = cam_out%dstdry3(i) + dst_fluxes(3)
320 0 : cam_out%dstdry4(i) = cam_out%dstdry4(i) + dst_fluxes(4)
321 :
322 : ! in rare cases, integrated deposition tendency is upward
323 0 : if (cam_out%bcphidry(i) < 0._r8) cam_out%bcphidry(i) = 0._r8
324 0 : if (cam_out%ocphidry(i) < 0._r8) cam_out%ocphidry(i) = 0._r8
325 0 : if (cam_out%bcphodry(i) < 0._r8) cam_out%bcphodry(i) = 0._r8
326 0 : if (cam_out%ocphodry(i) < 0._r8) cam_out%ocphodry(i) = 0._r8
327 0 : if (cam_out%dstdry1(i) < 0._r8) cam_out%dstdry1(i) = 0._r8
328 0 : if (cam_out%dstdry2(i) < 0._r8) cam_out%dstdry2(i) = 0._r8
329 0 : if (cam_out%dstdry3(i) < 0._r8) cam_out%dstdry3(i) = 0._r8
330 0 : if (cam_out%dstdry4(i) < 0._r8) cam_out%dstdry4(i) = 0._r8
331 :
332 : enddo
333 :
334 0 : end subroutine aero_deposition_cam_setdry
335 :
336 : end module aero_deposition_cam
|