Line data Source code
1 : module carma_aerosol_state_mod
2 : use shr_kind_mod, only: r8 => shr_kind_r8
3 : use aerosol_state_mod, only: aerosol_state, ptr2d_t
4 :
5 : use rad_constituents, only: rad_cnst_get_bin_mmr_by_idx, rad_cnst_get_bin_num !, rad_cnst_get_bin_mmr
6 : use rad_constituents, only: rad_cnst_get_info_by_bin
7 : use physics_buffer, only: physics_buffer_desc, pbuf_get_field, pbuf_get_index
8 : use physics_types, only: physics_state
9 : use aerosol_properties_mod, only: aerosol_properties, aero_name_len
10 :
11 : use physconst, only: pi
12 : use carma_intr, only: carma_get_total_mmr, carma_get_dry_radius, carma_get_number, carma_get_number_cld
13 : use carma_intr, only: carma_get_group_by_name, carma_get_kappa, carma_get_dry_radius, carma_get_wet_radius
14 : use carma_intr, only: carma_get_wght_pct
15 : use ppgrid, only: begchunk, endchunk, pcols, pver
16 :
17 : implicit none
18 :
19 : private
20 :
21 : public :: carma_aerosol_state
22 :
23 : type, extends(aerosol_state) :: carma_aerosol_state
24 : private
25 : type(physics_state), pointer :: state => null()
26 : type(physics_buffer_desc), pointer :: pbuf(:) => null()
27 : contains
28 :
29 : procedure :: get_transported
30 : procedure :: set_transported
31 : procedure :: ambient_total_bin_mmr
32 : procedure :: get_ambient_mmr_0list
33 : procedure :: get_ambient_mmr_rlist
34 : procedure :: get_cldbrne_mmr
35 : procedure :: get_ambient_num
36 : procedure :: get_cldbrne_num
37 : procedure :: get_states
38 : procedure :: icenuc_size_wght_arr
39 : procedure :: icenuc_size_wght_val
40 : procedure :: update_bin
41 : procedure :: hetfrz_size_wght
42 : procedure :: hygroscopicity
43 : procedure :: water_uptake
44 : procedure :: wgtpct
45 : procedure :: dry_volume
46 : procedure :: wet_volume
47 : procedure :: water_volume
48 : procedure :: wet_diameter
49 :
50 : final :: destructor
51 :
52 : end type carma_aerosol_state
53 :
54 : interface carma_aerosol_state
55 : procedure :: constructor
56 : end interface carma_aerosol_state
57 :
58 : real(r8), parameter :: four_thirds_pi = pi * 4._r8 / 3._r8
59 :
60 : contains
61 :
62 : !------------------------------------------------------------------------------
63 : !------------------------------------------------------------------------------
64 0 : function constructor(state,pbuf) result(newobj)
65 : type(physics_state), target, optional :: state
66 : type(physics_buffer_desc), pointer, optional :: pbuf(:)
67 :
68 : type(carma_aerosol_state), pointer :: newobj
69 :
70 : integer :: ierr
71 :
72 0 : allocate(newobj,stat=ierr)
73 0 : if( ierr /= 0 ) then
74 0 : nullify(newobj)
75 : return
76 : end if
77 :
78 0 : newobj%state => state
79 0 : newobj%pbuf => pbuf
80 :
81 0 : end function constructor
82 :
83 : !------------------------------------------------------------------------------
84 : !------------------------------------------------------------------------------
85 0 : subroutine destructor(self)
86 : type(carma_aerosol_state), intent(inout) :: self
87 :
88 0 : nullify(self%state)
89 0 : nullify(self%pbuf)
90 :
91 0 : end subroutine destructor
92 :
93 : !------------------------------------------------------------------------------
94 : ! sets transported components
95 : ! This aerosol model with the state of the transported aerosol constituents
96 : ! (mass mixing ratios or number mixing ratios)
97 : !------------------------------------------------------------------------------
98 0 : subroutine set_transported( self, transported_array )
99 : class(carma_aerosol_state), intent(inout) :: self
100 : real(r8), intent(in) :: transported_array(:,:,:)
101 : ! to be implemented later
102 0 : end subroutine set_transported
103 :
104 : !------------------------------------------------------------------------------
105 : ! returns transported components
106 : ! This returns to current state of the transported aerosol constituents
107 : ! (mass mixing ratios or number mixing ratios)
108 : !------------------------------------------------------------------------------
109 0 : subroutine get_transported( self, transported_array )
110 : class(carma_aerosol_state), intent(in) :: self
111 : real(r8), intent(out) :: transported_array(:,:,:)
112 : ! to be implemented later
113 0 : end subroutine get_transported
114 :
115 : !------------------------------------------------------------------------
116 : ! Total aerosol mass mixing ratio for a bin in a given grid box location (column and layer)
117 : !------------------------------------------------------------------------
118 0 : function ambient_total_bin_mmr(self, aero_props, bin_ndx, col_ndx, lyr_ndx) result(mmr_tot)
119 : class(carma_aerosol_state), intent(in) :: self
120 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
121 : integer, intent(in) :: bin_ndx ! bin index
122 : integer, intent(in) :: col_ndx ! column index
123 : integer, intent(in) :: lyr_ndx ! vertical layer index
124 :
125 : real(r8) :: mmr_tot ! mass mixing ratios totaled for all species
126 :
127 : real(r8) :: totmmr(pcols,pver)
128 : character(len=aero_name_len) :: bin_name, shortname
129 : integer :: igroup, ibin, rc, nchr
130 :
131 0 : call rad_cnst_get_info_by_bin(0, bin_ndx, bin_name=bin_name)
132 :
133 0 : nchr = len_trim(bin_name)-2
134 0 : shortname = bin_name(:nchr)
135 :
136 0 : call carma_get_group_by_name(shortname, igroup, rc)
137 :
138 0 : read(bin_name(nchr+1:),*) ibin
139 :
140 0 : call carma_get_total_mmr(self%state, igroup, ibin, totmmr, rc)
141 :
142 0 : mmr_tot = totmmr(col_ndx,lyr_ndx)
143 :
144 0 : end function ambient_total_bin_mmr
145 :
146 : !------------------------------------------------------------------------------
147 : ! returns ambient aerosol mass mixing ratio for a given species index and bin index
148 : !------------------------------------------------------------------------------
149 0 : subroutine get_ambient_mmr_0list(self, species_ndx, bin_ndx, mmr)
150 : class(carma_aerosol_state), intent(in) :: self
151 : integer, intent(in) :: species_ndx ! species index
152 : integer, intent(in) :: bin_ndx ! bin index
153 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
154 :
155 0 : call rad_cnst_get_bin_mmr_by_idx(0, bin_ndx, species_ndx, 'a', self%state, self%pbuf, mmr)
156 :
157 0 : end subroutine get_ambient_mmr_0list
158 :
159 : !------------------------------------------------------------------------------
160 : ! returns ambient aerosol mass mixing ratio for a given radiation diagnostics
161 : ! list index, species index and bin index
162 : !------------------------------------------------------------------------------
163 0 : subroutine get_ambient_mmr_rlist(self, list_ndx, species_ndx, bin_ndx, mmr)
164 : class(carma_aerosol_state), intent(in) :: self
165 : integer, intent(in) :: list_ndx ! rad climate list index
166 : integer, intent(in) :: species_ndx ! species index
167 : integer, intent(in) :: bin_ndx ! bin index
168 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
169 :
170 0 : call rad_cnst_get_bin_mmr_by_idx(list_ndx, bin_ndx, species_ndx, 'a', self%state, self%pbuf, mmr)
171 :
172 0 : end subroutine get_ambient_mmr_rlist
173 :
174 : !------------------------------------------------------------------------------
175 : ! returns cloud-borne aerosol number mixing ratio for a given species index and bin index
176 : !------------------------------------------------------------------------------
177 0 : subroutine get_cldbrne_mmr(self, species_ndx, bin_ndx, mmr)
178 : class(carma_aerosol_state), intent(in) :: self
179 : integer, intent(in) :: species_ndx ! species index
180 : integer, intent(in) :: bin_ndx ! bin index
181 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
182 :
183 0 : call rad_cnst_get_bin_mmr_by_idx(0, bin_ndx, species_ndx, 'c', self%state, self%pbuf, mmr)
184 :
185 0 : end subroutine get_cldbrne_mmr
186 :
187 : !------------------------------------------------------------------------------
188 : ! returns ambient aerosol number mixing ratio for a given species index and bin index
189 : !------------------------------------------------------------------------------
190 0 : subroutine get_ambient_num(self, bin_ndx, num)
191 : class(carma_aerosol_state), intent(in) :: self
192 : integer, intent(in) :: bin_ndx ! bin index
193 : real(r8), pointer :: num(:,:) ! number mixing ratios
194 :
195 : character(len=aero_name_len) :: bin_name, shortname
196 : integer :: igroup, ibin, rc, nchr, ncol
197 : real(r8) :: nmr(pcols,pver)
198 :
199 0 : ncol = self%state%ncol
200 :
201 0 : call rad_cnst_get_info_by_bin(0, bin_ndx, bin_name=bin_name)
202 :
203 0 : nchr = len_trim(bin_name)-2
204 0 : shortname = bin_name(:nchr)
205 :
206 0 : call carma_get_group_by_name(shortname, igroup, rc)
207 :
208 0 : read(bin_name(nchr+1:),*) ibin
209 :
210 0 : call rad_cnst_get_bin_num(0, bin_ndx, 'a', self%state, self%pbuf, num)
211 :
212 0 : call carma_get_number(self%state, igroup, ibin, nmr, rc)
213 :
214 0 : num(:ncol,:) = nmr(:ncol,:)
215 :
216 0 : end subroutine get_ambient_num
217 :
218 : !------------------------------------------------------------------------------
219 : ! returns cloud-borne aerosol number mixing ratio for a given species index and bin index
220 : !------------------------------------------------------------------------------
221 0 : subroutine get_cldbrne_num(self, bin_ndx, num)
222 : class(carma_aerosol_state), intent(in) :: self
223 : integer, intent(in) :: bin_ndx ! bin index
224 : real(r8), pointer :: num(:,:) ! number mixing ratios
225 :
226 : character(len=aero_name_len) :: bin_name, shortname
227 : integer :: igroup, ibin, rc, nchr, ncol
228 : real(r8) :: nmr(pcols,pver)
229 :
230 0 : ncol = self%state%ncol
231 :
232 0 : call rad_cnst_get_info_by_bin(0, bin_ndx, bin_name=bin_name)
233 :
234 0 : nchr = len_trim(bin_name)-2
235 0 : shortname = bin_name(:nchr)
236 :
237 0 : call carma_get_group_by_name(shortname, igroup, rc)
238 :
239 0 : read(bin_name(nchr+1:),*) ibin
240 :
241 0 : call rad_cnst_get_bin_num(0, bin_ndx, 'c', self%state, self%pbuf, num)
242 :
243 0 : call carma_get_number_cld(self%pbuf, igroup, ibin, ncol, pver, nmr, rc)
244 :
245 0 : num(:ncol,:) = nmr(:ncol,:)
246 :
247 0 : end subroutine get_cldbrne_num
248 :
249 : !------------------------------------------------------------------------------
250 : ! returns interstitial and cloud-borne aerosol states
251 : !------------------------------------------------------------------------------
252 0 : subroutine get_states( self, aero_props, raer, qqcw )
253 : class(carma_aerosol_state), intent(in) :: self
254 : class(aerosol_properties), intent(in) :: aero_props
255 : type(ptr2d_t), intent(out) :: raer(:)
256 : type(ptr2d_t), intent(out) :: qqcw(:)
257 :
258 : integer :: ibin,ispc, indx
259 :
260 0 : do ibin = 1, aero_props%nbins()
261 0 : indx = aero_props%indexer(ibin, 0)
262 0 : call self%get_ambient_num(ibin, raer(indx)%fld)
263 0 : call self%get_cldbrne_num(ibin, qqcw(indx)%fld)
264 0 : do ispc = 1, aero_props%nspecies(ibin)
265 0 : indx = aero_props%indexer(ibin, ispc)
266 0 : call self%get_ambient_mmr(ispc,ibin, raer(indx)%fld)
267 0 : call self%get_cldbrne_mmr(ispc,ibin, qqcw(indx)%fld)
268 : end do
269 : end do
270 :
271 0 : end subroutine get_states
272 :
273 : !------------------------------------------------------------------------------
274 : ! return aerosol bin size weights for a given bin
275 : !------------------------------------------------------------------------------
276 0 : subroutine icenuc_size_wght_arr(self, bin_ndx, ncol, nlev, species_type, use_preexisting_ice, wght)
277 : class(carma_aerosol_state), intent(in) :: self
278 : integer, intent(in) :: bin_ndx ! bin number
279 : integer, intent(in) :: ncol ! number of columns
280 : integer, intent(in) :: nlev ! number of vertical levels
281 : character(len=*), intent(in) :: species_type ! species type
282 : logical, intent(in) :: use_preexisting_ice ! pre-existing ice flag
283 : real(r8), intent(out) :: wght(:,:)
284 :
285 : character(len=aero_name_len) :: bin_name, shortname
286 0 : real(r8) :: rdry(ncol,nlev), rhopdry(ncol,nlev)
287 : integer :: i,k
288 : real(r8) :: diamdry
289 : integer :: igroup, ibin, rc, nchr
290 :
291 0 : wght = 0._r8
292 :
293 0 : call rad_cnst_get_info_by_bin(0, bin_ndx, bin_name=bin_name)
294 :
295 0 : nchr = len_trim(bin_name)-2
296 0 : shortname = bin_name(:nchr)
297 0 : call carma_get_group_by_name(shortname, igroup, rc)
298 :
299 0 : read(bin_name(nchr+1:),*) ibin
300 :
301 0 : call carma_get_dry_radius(self%state, igroup, ibin, rdry, rhopdry, rc) ! m, kg/m3
302 :
303 0 : do k = 1,nlev
304 0 : do i = 1,ncol
305 0 : diamdry = rdry(i,k) * 2._r8 * 1.e6_r8 ! diameter in microns (from radius in m)
306 0 : if (diamdry >= 0.1_r8) then ! size threashold
307 0 : wght(i,k) = 1._r8
308 : end if
309 : end do
310 : end do
311 :
312 0 : end subroutine icenuc_size_wght_arr
313 :
314 : !------------------------------------------------------------------------------
315 : ! return aerosol bin size weights for a given bin, column and vertical layer
316 : !------------------------------------------------------------------------------
317 0 : subroutine icenuc_size_wght_val(self, bin_ndx, col_ndx, lyr_ndx, species_type, use_preexisting_ice, wght)
318 : class(carma_aerosol_state), intent(in) :: self
319 : integer, intent(in) :: bin_ndx ! bin number
320 : integer, intent(in) :: col_ndx ! column index
321 : integer, intent(in) :: lyr_ndx ! vertical layer index
322 : character(len=*), intent(in) :: species_type ! species type
323 : logical, intent(in) :: use_preexisting_ice ! pre-existing ice flag
324 : real(r8), intent(out) :: wght
325 :
326 : real(r8) :: wght_arr(pcols,pver)
327 :
328 0 : call self%icenuc_size_wght(bin_ndx, self%state%ncol, pver, species_type, use_preexisting_ice, wght_arr)
329 :
330 0 : wght = wght_arr(col_ndx,lyr_ndx)
331 :
332 0 : end subroutine icenuc_size_wght_val
333 :
334 : !------------------------------------------------------------------------------
335 : !------------------------------------------------------------------------------
336 0 : subroutine update_bin( self, bin_ndx, col_ndx, lyr_ndx, delmmr_sum, delnum_sum, tnd_ndx, dtime, tend )
337 : class(carma_aerosol_state), intent(in) :: self
338 : integer, intent(in) :: bin_ndx ! bin number
339 : integer, intent(in) :: col_ndx ! column index
340 : integer, intent(in) :: lyr_ndx ! vertical layer index
341 : real(r8),intent(in) :: delmmr_sum ! mass mixing ratio change summed over all species in bin
342 : real(r8),intent(in) :: delnum_sum ! number mixing ratio change summed over all species in bin
343 : integer, intent(in) :: tnd_ndx ! tendency index
344 : real(r8),intent(in) :: dtime ! time step size (sec)
345 : real(r8),intent(inout) :: tend(:,:,:) ! tendency
346 :
347 : real(r8), pointer :: amb_num(:,:)
348 : real(r8), pointer :: cld_num(:,:)
349 :
350 : ! for updating num (num tendancies)
351 : ! -- nothing to do here for CARMA since num is calculated when needed
352 :
353 0 : end subroutine update_bin
354 :
355 : !------------------------------------------------------------------------------
356 : ! returns the volume-weighted fractions of aerosol subset `bin_ndx` that can act
357 : ! as heterogeneous freezing nuclei
358 : !------------------------------------------------------------------------------
359 0 : function hetfrz_size_wght(self, bin_ndx, ncol, nlev) result(wght)
360 : class(carma_aerosol_state), intent(in) :: self
361 : integer, intent(in) :: bin_ndx ! bin number
362 : integer, intent(in) :: ncol ! number of columns
363 : integer, intent(in) :: nlev ! number of vertical levels
364 :
365 : real(r8) :: wght(ncol,nlev)
366 :
367 : character(len=aero_name_len) :: bin_name, shortname
368 0 : real(r8) :: rdry(ncol,nlev), rhopdry(ncol,nlev)
369 : integer :: i,k
370 : real(r8) :: diamdry
371 : integer :: igroup, ibin, rc, nchr
372 :
373 0 : wght = 0._r8
374 :
375 0 : call rad_cnst_get_info_by_bin(0, bin_ndx, bin_name=bin_name)
376 :
377 0 : nchr = len_trim(bin_name)-2
378 0 : shortname = bin_name(:nchr)
379 0 : call carma_get_group_by_name(shortname, igroup, rc)
380 :
381 0 : read(bin_name(nchr+1:),*) ibin
382 :
383 0 : call carma_get_dry_radius(self%state, igroup, ibin, rdry, rhopdry, rc) ! m, kg/m3
384 :
385 0 : do k = 1,nlev
386 0 : do i = 1,ncol
387 0 : diamdry = rdry(i,k) * 2._r8 * 1.e6_r8 ! diameter in microns (from radius in m)
388 0 : if (diamdry >= 0.1_r8) then ! size threashold
389 0 : wght(i,k) = 1._r8
390 : end if
391 : end do
392 : end do
393 :
394 0 : end function hetfrz_size_wght
395 :
396 : !------------------------------------------------------------------------------
397 : ! returns hygroscopicity for a given radiation diagnostic list number and
398 : ! bin number
399 : !------------------------------------------------------------------------------
400 0 : subroutine hygroscopicity(self, list_ndx, bin_ndx, kappa)
401 : class(carma_aerosol_state), intent(in) :: self
402 : integer, intent(in) :: list_ndx ! rad climate list number
403 : integer, intent(in) :: bin_ndx ! bin number
404 : real(r8), intent(out) :: kappa(:,:) ! hygroscopicity (ncol,nlev)
405 :
406 : character(len=aero_name_len) :: bin_name, shortname
407 : integer :: igroup, ibin, rc, nchr, ncol
408 :
409 0 : call rad_cnst_get_info_by_bin(0, bin_ndx, bin_name=bin_name)
410 :
411 0 : nchr = len_trim(bin_name)-2
412 0 : shortname = bin_name(:nchr)
413 :
414 0 : call carma_get_group_by_name(shortname, igroup, rc)
415 :
416 0 : read(bin_name(nchr+1:),*) ibin
417 :
418 0 : call carma_get_kappa(self%state, igroup, ibin, kappa, rc)
419 :
420 0 : end subroutine hygroscopicity
421 :
422 : !------------------------------------------------------------------------------
423 : ! returns aerosol wet diameter and aerosol water concentration for a given
424 : ! radiation diagnostic list number and bin number
425 : !------------------------------------------------------------------------------
426 0 : subroutine water_uptake(self, aero_props, list_idx, bin_idx, ncol, nlev, dgnumwet, qaerwat)
427 :
428 : class(carma_aerosol_state), intent(in) :: self
429 : class(aerosol_properties), intent(in) :: aero_props
430 : integer, intent(in) :: list_idx ! rad climate/diags list number
431 : integer, intent(in) :: bin_idx ! bin number
432 : integer, intent(in) :: ncol ! number of columns
433 : integer, intent(in) :: nlev ! number of levels
434 : real(r8),intent(out) :: dgnumwet(ncol,nlev) ! aerosol wet diameter (m)
435 : real(r8),intent(out) :: qaerwat(ncol,nlev) ! aerosol water concentration (g/g)
436 :
437 0 : dgnumwet = -huge(1._r8)
438 0 : qaerwat = -huge(1._r8)
439 :
440 0 : end subroutine water_uptake
441 :
442 : !------------------------------------------------------------------------------
443 : ! aerosol weight precent of H2SO4/H2O solution
444 : !------------------------------------------------------------------------------
445 0 : function wgtpct(self, ncol, nlev) result(wtp)
446 : class(carma_aerosol_state), intent(in) :: self
447 : integer, intent(in) :: ncol, nlev
448 : real(r8) :: wtp(ncol,nlev) ! weight precent of H2SO4/H2O solution for given icol, ilev
449 :
450 0 : wtp(:,:) = carma_get_wght_pct(ncol,nlev,self%state)
451 :
452 0 : end function wgtpct
453 :
454 : !------------------------------------------------------------------------------
455 : ! aerosol dry volume (m3/kg) for given radiation diagnostic list number and bin number
456 : !------------------------------------------------------------------------------
457 0 : function dry_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
458 :
459 : class(carma_aerosol_state), intent(in) :: self
460 : class(aerosol_properties), intent(in) :: aero_props
461 :
462 : integer, intent(in) :: list_idx ! rad climate/diags list number
463 : integer, intent(in) :: bin_idx ! bin number
464 : integer, intent(in) :: ncol ! number of columns
465 : integer, intent(in) :: nlev ! number of levels
466 :
467 : real(r8) :: vol(ncol,nlev) ! m3/kg
468 :
469 : real(r8) :: raddry(pcols,pver) ! dry radius (m)
470 : real(r8) :: rhodry(pcols,pver) ! dry density (kg/m3)
471 : real(r8) :: nmr(pcols,pver) ! number mixing ratio (#/kg)
472 :
473 : character(len=aero_name_len) :: bin_name, shortname
474 : integer :: igroup, ibin, rc, nchr
475 :
476 0 : call rad_cnst_get_info_by_bin(0, bin_idx, bin_name=bin_name)
477 :
478 0 : nchr = len_trim(bin_name)-2
479 0 : shortname = bin_name(:nchr)
480 :
481 0 : call carma_get_group_by_name(shortname, igroup, rc)
482 :
483 0 : read(bin_name(nchr+1:),*) ibin
484 :
485 0 : vol = 0._r8
486 :
487 0 : call carma_get_dry_radius(self%state, igroup, ibin, raddry, rhodry, rc)
488 0 : call carma_get_number(self%state, igroup, ibin, nmr, rc)
489 :
490 0 : vol(:ncol,:) = four_thirds_pi * (raddry(:ncol,:)**3) * nmr(:ncol,:) ! units = m3/kg
491 :
492 0 : end function dry_volume
493 :
494 : !------------------------------------------------------------------------------
495 : ! aerosol wet volume (m3/kg) for given radiation diagnostic list number and bin number
496 : !------------------------------------------------------------------------------
497 0 : function wet_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
498 :
499 : class(carma_aerosol_state), intent(in) :: self
500 : class(aerosol_properties), intent(in) :: aero_props
501 :
502 : integer, intent(in) :: list_idx ! rad climate/diags list number
503 : integer, intent(in) :: bin_idx ! bin number
504 : integer, intent(in) :: ncol ! number of columns
505 : integer, intent(in) :: nlev ! number of levels
506 :
507 : real(r8) :: vol(ncol,nlev) ! m3/kg
508 :
509 : real(r8) :: radwet(pcols,pver) ! wet radius (m)
510 : real(r8) :: rhowet(pcols,pver) ! wet density (kg/m3)
511 : real(r8) :: nmr(pcols,pver) ! number mixing ratio (#/kg)
512 :
513 : character(len=aero_name_len) :: bin_name, shortname
514 : integer :: igroup, ibin, rc, nchr
515 :
516 0 : call rad_cnst_get_info_by_bin(0, bin_idx, bin_name=bin_name)
517 :
518 0 : nchr = len_trim(bin_name)-2
519 0 : shortname = bin_name(:nchr)
520 :
521 0 : call carma_get_group_by_name(shortname, igroup, rc)
522 :
523 0 : read(bin_name(nchr+1:),*) ibin
524 :
525 0 : vol = 0._r8
526 :
527 0 : call carma_get_wet_radius(self%state, igroup, ibin, radwet, rhowet, rc)
528 0 : call carma_get_number(self%state, igroup, ibin, nmr, rc)
529 :
530 0 : vol(:ncol,:) = four_thirds_pi * (radwet(:ncol,:)**3) * nmr(:ncol,:) ! units = m3/kg
531 :
532 0 : end function wet_volume
533 :
534 : !------------------------------------------------------------------------------
535 : ! aerosol water volume (m3/kg) for given radiation diagnostic list number and bin number
536 : !------------------------------------------------------------------------------
537 0 : function water_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
538 :
539 : class(carma_aerosol_state), intent(in) :: self
540 : class(aerosol_properties), intent(in) :: aero_props
541 :
542 : integer, intent(in) :: list_idx ! rad climate/diags list number
543 : integer, intent(in) :: bin_idx ! bin number
544 : integer, intent(in) :: ncol ! number of columns
545 : integer, intent(in) :: nlev ! number of levels
546 :
547 : real(r8) :: vol(ncol,nlev) ! m3/kg
548 :
549 0 : real(r8) :: wetvol(ncol,nlev)
550 0 : real(r8) :: dryvol(ncol,nlev)
551 :
552 0 : wetvol = self%wet_volume(aero_props, list_idx, bin_idx, ncol, nlev)
553 0 : dryvol = self%dry_volume(aero_props, list_idx, bin_idx, ncol, nlev)
554 :
555 0 : vol(:ncol,:) = wetvol(:ncol,:) - dryvol(:ncol,:)
556 :
557 0 : where (vol<0._r8)
558 : vol = 0._r8
559 : end where
560 :
561 0 : end function water_volume
562 :
563 : !------------------------------------------------------------------------------
564 : ! aerosol wet diameter
565 : !------------------------------------------------------------------------------
566 0 : function wet_diameter(self, bin_idx, ncol, nlev) result(diam)
567 : class(carma_aerosol_state), intent(in) :: self
568 : integer, intent(in) :: bin_idx ! bin number
569 : integer, intent(in) :: ncol ! number of columns
570 : integer, intent(in) :: nlev ! number of levels
571 :
572 : real(r8) :: diam(ncol,nlev)
573 :
574 : real(r8) :: radwet(pcols,pver) !! wet radius (m)
575 : real(r8) :: rhowet(pcols,pver) !! wet density (kg/m3)
576 :
577 : character(len=aero_name_len) :: bin_name, shortname
578 : integer :: igroup, ibin, rc, nchr
579 :
580 0 : call rad_cnst_get_info_by_bin(0, bin_idx, bin_name=bin_name)
581 :
582 0 : nchr = len_trim(bin_name)-2
583 0 : shortname = bin_name(:nchr)
584 :
585 0 : call carma_get_group_by_name(shortname, igroup, rc)
586 :
587 0 : read(bin_name(nchr+1:),*) ibin
588 :
589 0 : call carma_get_wet_radius(self%state, igroup, ibin, radwet, rhowet, rc)
590 :
591 0 : diam(:ncol,:nlev) = 2._r8*radwet(:ncol,:nlev)
592 :
593 0 : end function wet_diameter
594 :
595 0 : end module carma_aerosol_state_mod
|