Line data Source code
1 : module modal_aerosol_state_mod
2 : use shr_kind_mod, only: r8 => shr_kind_r8
3 : use shr_spfn_mod, only: erf => shr_spfn_erf
4 : use aerosol_state_mod, only: aerosol_state, ptr2d_t
5 : use rad_constituents, only: rad_cnst_get_aer_mmr, rad_cnst_get_mode_num, rad_cnst_get_info
6 : use rad_constituents, only: rad_cnst_get_mode_props
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 : use physconst, only: rhoh2o
11 :
12 : implicit none
13 :
14 : private
15 :
16 : public :: modal_aerosol_state
17 :
18 : type, extends(aerosol_state) :: modal_aerosol_state
19 : private
20 : type(physics_state), pointer :: state => null()
21 : type(physics_buffer_desc), pointer :: pbuf(:) => null()
22 : contains
23 :
24 : procedure :: get_transported
25 : procedure :: set_transported
26 : procedure :: ambient_total_bin_mmr
27 : procedure :: get_ambient_mmr_0list
28 : procedure :: get_ambient_mmr_rlist
29 : procedure :: get_cldbrne_mmr
30 : procedure :: get_ambient_num
31 : procedure :: get_cldbrne_num
32 : procedure :: get_states
33 : procedure :: icenuc_size_wght_arr
34 : procedure :: icenuc_size_wght_val
35 : procedure :: icenuc_type_wght
36 : procedure :: update_bin
37 : procedure :: hetfrz_size_wght
38 : procedure :: hygroscopicity
39 : procedure :: water_uptake
40 : procedure :: dry_volume
41 : procedure :: wet_volume
42 : procedure :: water_volume
43 :
44 : final :: destructor
45 :
46 : end type modal_aerosol_state
47 :
48 : interface modal_aerosol_state
49 : procedure :: constructor
50 : end interface modal_aerosol_state
51 :
52 : real(r8), parameter :: rh2odens = 1._r8/rhoh2o
53 :
54 : contains
55 :
56 : !------------------------------------------------------------------------------
57 : !------------------------------------------------------------------------------
58 244584 : function constructor(state,pbuf) result(newobj)
59 : type(physics_state), target :: state
60 : type(physics_buffer_desc), pointer :: pbuf(:)
61 :
62 : type(modal_aerosol_state), pointer :: newobj
63 :
64 : integer :: ierr
65 :
66 244584 : allocate(newobj,stat=ierr)
67 244584 : if( ierr /= 0 ) then
68 244584 : nullify(newobj)
69 : return
70 : end if
71 :
72 244584 : newobj%state => state
73 244584 : newobj%pbuf => pbuf
74 :
75 244584 : end function constructor
76 :
77 : !------------------------------------------------------------------------------
78 : !------------------------------------------------------------------------------
79 244584 : subroutine destructor(self)
80 : type(modal_aerosol_state), intent(inout) :: self
81 :
82 244584 : nullify(self%state)
83 244584 : nullify(self%pbuf)
84 :
85 244584 : end subroutine destructor
86 :
87 : !------------------------------------------------------------------------------
88 : ! sets transported components
89 : ! This aerosol model with the state of the transported aerosol constituents
90 : ! (mass mixing ratios or number mixing ratios)
91 : !------------------------------------------------------------------------------
92 0 : subroutine set_transported( self, transported_array )
93 : class(modal_aerosol_state), intent(inout) :: self
94 : real(r8), intent(in) :: transported_array(:,:,:)
95 : ! to be implemented later
96 0 : end subroutine set_transported
97 :
98 : !------------------------------------------------------------------------------
99 : ! returns transported components
100 : ! This returns to current state of the transported aerosol constituents
101 : ! (mass mixing ratios or number mixing ratios)
102 : !------------------------------------------------------------------------------
103 0 : subroutine get_transported( self, transported_array )
104 : class(modal_aerosol_state), intent(in) :: self
105 : real(r8), intent(out) :: transported_array(:,:,:)
106 : ! to be implemented later
107 0 : end subroutine get_transported
108 :
109 : !------------------------------------------------------------------------
110 : ! Total aerosol mass mixing ratio for a bin in a given grid box location (column and layer)
111 : !------------------------------------------------------------------------
112 0 : function ambient_total_bin_mmr(self, aero_props, bin_ndx, col_ndx, lyr_ndx) result(mmr_tot)
113 : class(modal_aerosol_state), intent(in) :: self
114 : class(aerosol_properties), intent(in) :: aero_props
115 : integer, intent(in) :: bin_ndx ! bin index
116 : integer, intent(in) :: col_ndx ! column index
117 : integer, intent(in) :: lyr_ndx ! vertical layer index
118 :
119 : real(r8) :: mmr_tot ! mass mixing ratios totaled for all species
120 0 : real(r8),pointer :: mmrptr(:,:)
121 : integer :: spec_ndx
122 :
123 0 : mmr_tot = 0._r8
124 :
125 0 : do spec_ndx=1,aero_props%nspecies(bin_ndx)
126 0 : call rad_cnst_get_aer_mmr(0, bin_ndx, spec_ndx, 'a', self%state, self%pbuf, mmrptr)
127 0 : mmr_tot = mmr_tot + mmrptr(col_ndx,lyr_ndx)
128 : end do
129 :
130 0 : end function ambient_total_bin_mmr
131 :
132 : !------------------------------------------------------------------------------
133 : ! returns ambient aerosol mass mixing ratio for a given species index and bin index
134 : !------------------------------------------------------------------------------
135 627727017 : subroutine get_ambient_mmr_0list(self, species_ndx, bin_ndx, mmr)
136 : class(modal_aerosol_state), intent(in) :: self
137 : integer, intent(in) :: species_ndx ! species index
138 : integer, intent(in) :: bin_ndx ! bin index
139 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
140 :
141 627727017 : call rad_cnst_get_aer_mmr(0, bin_ndx, species_ndx, 'a', self%state, self%pbuf, mmr)
142 627727017 : end subroutine get_ambient_mmr_0list
143 :
144 : !------------------------------------------------------------------------------
145 : ! returns ambient aerosol mass mixing ratio for a given radiation diagnostics
146 : ! list index, species index and bin index
147 : !------------------------------------------------------------------------------
148 2060488800 : subroutine get_ambient_mmr_rlist(self, list_ndx, species_ndx, bin_ndx, mmr)
149 : class(modal_aerosol_state), intent(in) :: self
150 : integer, intent(in) :: list_ndx ! rad climate list index
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 2060488800 : call rad_cnst_get_aer_mmr(list_ndx, bin_ndx, species_ndx, 'a', self%state, self%pbuf, mmr)
156 2060488800 : end subroutine get_ambient_mmr_rlist
157 :
158 : !------------------------------------------------------------------------------
159 : ! returns cloud-borne aerosol number mixing ratio for a given species index and bin index
160 : !------------------------------------------------------------------------------
161 614491617 : subroutine get_cldbrne_mmr(self, species_ndx, bin_ndx, mmr)
162 : class(modal_aerosol_state), intent(in) :: self
163 : integer, intent(in) :: species_ndx ! species index
164 : integer, intent(in) :: bin_ndx ! bin index
165 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
166 :
167 614491617 : call rad_cnst_get_aer_mmr(0, bin_ndx, species_ndx, 'c', self%state, self%pbuf, mmr)
168 614491617 : end subroutine get_cldbrne_mmr
169 :
170 : !------------------------------------------------------------------------------
171 : ! returns ambient aerosol number mixing ratio for a given species index and bin index
172 : !------------------------------------------------------------------------------
173 479377172 : subroutine get_ambient_num(self, bin_ndx, num)
174 : class(modal_aerosol_state), intent(in) :: self
175 : integer, intent(in) :: bin_ndx ! bin index
176 : real(r8), pointer :: num(:,:) ! number densities
177 :
178 479377172 : call rad_cnst_get_mode_num(0, bin_ndx, 'a', self%state, self%pbuf, num)
179 479377172 : end subroutine get_ambient_num
180 :
181 : !------------------------------------------------------------------------------
182 : ! returns cloud-borne aerosol number mixing ratio for a given species index and bin index
183 : !------------------------------------------------------------------------------
184 476553620 : subroutine get_cldbrne_num(self, bin_ndx, num)
185 : class(modal_aerosol_state), intent(in) :: self
186 : integer, intent(in) :: bin_ndx ! bin index
187 : real(r8), pointer :: num(:,:)
188 :
189 476553620 : call rad_cnst_get_mode_num(0, bin_ndx, 'c', self%state, self%pbuf, num)
190 476553620 : end subroutine get_cldbrne_num
191 :
192 : !------------------------------------------------------------------------------
193 : ! returns interstitial and cloud-borne aerosol states
194 : !------------------------------------------------------------------------------
195 176472 : subroutine get_states( self, aero_props, raer, qqcw )
196 : class(modal_aerosol_state), intent(in) :: self
197 : class(aerosol_properties), intent(in) :: aero_props
198 : type(ptr2d_t), intent(out) :: raer(:)
199 : type(ptr2d_t), intent(out) :: qqcw(:)
200 :
201 : integer :: ibin,ispc, indx
202 :
203 882360 : do ibin = 1, aero_props%nbins()
204 705888 : indx = aero_props%indexer(ibin, 0)
205 705888 : call self%get_ambient_num(ibin, raer(indx)%fld)
206 705888 : call self%get_cldbrne_num(ibin, qqcw(indx)%fld)
207 4235328 : do ispc = 1, aero_props%nspecies(ibin)
208 2647080 : indx = aero_props%indexer(ibin, ispc)
209 2647080 : call self%get_ambient_mmr(ispc,ibin, raer(indx)%fld)
210 3352968 : call self%get_cldbrne_mmr(ispc,ibin, qqcw(indx)%fld)
211 : end do
212 : end do
213 :
214 176472 : end subroutine get_states
215 :
216 : !------------------------------------------------------------------------------
217 : ! return aerosol bin size weights for a given bin
218 : !------------------------------------------------------------------------------
219 3352968 : subroutine icenuc_size_wght_arr(self, bin_ndx, ncol, nlev, species_type, use_preexisting_ice, wght)
220 : class(modal_aerosol_state), intent(in) :: self
221 : integer, intent(in) :: bin_ndx ! bin number
222 : integer, intent(in) :: ncol ! number of columns
223 : integer, intent(in) :: nlev ! number of vertical levels
224 : character(len=*), intent(in) :: species_type ! species type
225 : logical, intent(in) :: use_preexisting_ice ! pre-existing ice flag
226 : real(r8), intent(out) :: wght(:,:)
227 :
228 : character(len=aero_name_len) :: modetype
229 3352968 : real(r8), pointer :: dgnum(:,:,:) ! mode dry radius
230 : real(r8) :: sigmag_aitken
231 : integer :: i,k
232 :
233 3352968 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
234 :
235 5210122392 : wght = 0._r8
236 :
237 6705936 : select case ( trim(species_type) )
238 : case('dust')
239 529416 : if (modetype=='coarse' .or. modetype=='coarse_dust') then
240 274216968 : wght(:ncol,:) = 1._r8
241 : end if
242 : case('sulfate')
243 529416 : if (modetype=='aitken') then
244 176472 : if ( use_preexisting_ice ) then
245 274216968 : wght(:ncol,:) = 1._r8
246 : else
247 0 : call rad_cnst_get_mode_props(0, bin_ndx, sigmag=sigmag_aitken)
248 0 : call pbuf_get_field(self%pbuf, pbuf_get_index('DGNUM' ), dgnum)
249 0 : do k = 1,nlev
250 0 : do i = 1,ncol
251 0 : if (dgnum(i,k,bin_ndx) > 0._r8) then
252 : ! only allow so4 with D>0.1 um in ice nucleation
253 0 : wght(i,k) = max(0._r8,(0.5_r8 - 0.5_r8* &
254 : erf(log(0.1e-6_r8/dgnum(i,k,bin_ndx))/ &
255 0 : (2._r8**0.5_r8*log(sigmag_aitken))) ))
256 : end if
257 : end do
258 : end do
259 : endif
260 : endif
261 : case('black-c')
262 352944 : if (modetype=='accum') then
263 274216968 : wght(:ncol,:) = 1._r8
264 : endif
265 : case('sulfate_strat')
266 6705936 : if (modetype=='accum' .or. modetype=='coarse' .or. modetype=='coarse_strat') then
267 548433936 : wght(:ncol,:) = 1._r8
268 : endif
269 : end select
270 :
271 3352968 : end subroutine icenuc_size_wght_arr
272 :
273 : !------------------------------------------------------------------------------
274 : ! return aerosol bin size weights for a given bin, column and vertical layer
275 : !------------------------------------------------------------------------------
276 180451908 : subroutine icenuc_size_wght_val(self, bin_ndx, col_ndx, lyr_ndx, species_type, use_preexisting_ice, wght)
277 : class(modal_aerosol_state), intent(in) :: self
278 : integer, intent(in) :: bin_ndx ! bin number
279 : integer, intent(in) :: col_ndx ! column index
280 : integer, intent(in) :: lyr_ndx ! vertical layer index
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) :: modetype
286 180451908 : real(r8), pointer :: dgnum(:,:,:) ! mode dry radius
287 : real(r8) :: sigmag_aitken
288 :
289 180451908 : wght = 0._r8
290 :
291 180451908 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
292 :
293 360903816 : select case ( trim(species_type) )
294 : case('dust')
295 180451908 : if (modetype=='coarse' .or. modetype=='coarse_dust') then
296 180451908 : wght = 1._r8
297 : end if
298 : case('sulfate')
299 0 : if (modetype=='aitken') then
300 0 : if ( use_preexisting_ice ) then
301 0 : wght = 1._r8
302 : else
303 0 : call rad_cnst_get_mode_props(0, bin_ndx, sigmag=sigmag_aitken)
304 0 : call pbuf_get_field(self%pbuf, pbuf_get_index('DGNUM' ), dgnum)
305 :
306 0 : if (dgnum(col_ndx,lyr_ndx,bin_ndx) > 0._r8) then
307 : ! only allow so4 with D>0.1 um in ice nucleation
308 : wght = max(0._r8,(0.5_r8 - 0.5_r8* &
309 : erf(log(0.1e-6_r8/dgnum(col_ndx,lyr_ndx,bin_ndx))/ &
310 0 : (2._r8**0.5_r8*log(sigmag_aitken))) ))
311 :
312 : end if
313 : endif
314 : endif
315 : case('black-c')
316 0 : if (modetype=='accum') then
317 0 : wght = 1._r8
318 : endif
319 : case('sulfate_strat')
320 360903816 : if (modetype=='accum' .or. modetype=='coarse' .or. modetype=='coarse_strat') then
321 0 : wght = 1._r8
322 : endif
323 : end select
324 :
325 180451908 : end subroutine icenuc_size_wght_val
326 :
327 : !------------------------------------------------------------------------------
328 : ! returns aerosol type weights for a given aerosol type and bin
329 : !------------------------------------------------------------------------------
330 3352968 : subroutine icenuc_type_wght(self, bin_ndx, ncol, nlev, species_type, aero_props, rho, wght, cloud_borne)
331 :
332 : use aerosol_properties_mod, only: aerosol_properties
333 :
334 : class(modal_aerosol_state), intent(in) :: self
335 : integer, intent(in) :: bin_ndx ! bin number
336 : integer, intent(in) :: ncol ! number of columns
337 : integer, intent(in) :: nlev ! number of vertical levels
338 : character(len=*), intent(in) :: species_type ! species type
339 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
340 : real(r8), intent(in) :: rho(:,:) ! air density (kg m-3)
341 : real(r8), intent(out) :: wght(:,:) ! type weights
342 : logical, optional, intent(in) :: cloud_borne ! if TRUE cloud-borne aerosols are used
343 : ! otherwise ambient aerosols are used
344 :
345 : character(len=aero_name_len) :: modetype
346 :
347 3352968 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
348 :
349 5210122392 : wght = 0._r8
350 :
351 3352968 : if (species_type == 'dust') then
352 529416 : if (modetype=='coarse_dust') then
353 0 : wght(:ncol,:) = 1._r8
354 : else
355 529416 : call self%icenuc_type_wght_base(bin_ndx, ncol, nlev, species_type, aero_props, rho, wght, cloud_borne)
356 : end if
357 2823552 : else if (species_type == 'sulfate_strat') then
358 705888 : if (modetype=='accum') then
359 274216968 : wght(:ncol,:) = 1._r8
360 529416 : elseif ( modetype=='coarse' .or. modetype=='coarse_strat') then
361 176472 : call self%icenuc_type_wght_base(bin_ndx, ncol, nlev, species_type, aero_props, rho, wght, cloud_borne)
362 : endif
363 : else
364 3290603616 : wght(:ncol,:) = 1._r8
365 : end if
366 :
367 3352968 : end subroutine icenuc_type_wght
368 :
369 : !------------------------------------------------------------------------------
370 : !------------------------------------------------------------------------------
371 180451908 : subroutine update_bin( self, bin_ndx, col_ndx, lyr_ndx, delmmr_sum, delnum_sum, tnd_ndx, dtime, tend )
372 : class(modal_aerosol_state), intent(in) :: self
373 : integer, intent(in) :: bin_ndx ! bin number
374 : integer, intent(in) :: col_ndx ! column index
375 : integer, intent(in) :: lyr_ndx ! vertical layer index
376 : real(r8),intent(in) :: delmmr_sum ! mass mixing ratio change summed over all species in bin
377 : real(r8),intent(in) :: delnum_sum ! number mixing ratio change summed over all species in bin
378 : integer, intent(in) :: tnd_ndx ! tendency index
379 : real(r8),intent(in) :: dtime ! time step size (sec)
380 : real(r8),intent(inout) :: tend(:,:,:) ! tendency
381 :
382 180451908 : real(r8), pointer :: amb_num(:,:)
383 180451908 : real(r8), pointer :: cld_num(:,:)
384 :
385 180451908 : call self%get_ambient_num(bin_ndx, amb_num)
386 180451908 : call self%get_cldbrne_num(bin_ndx, cld_num)
387 :
388 : ! if there is no bin mass compute updates/tendencies for bin number
389 : ! -- apply the total number change to bin number
390 180451908 : if (tnd_ndx>0) then
391 180451908 : tend(col_ndx,lyr_ndx,tnd_ndx) = -delnum_sum/dtime
392 : else
393 0 : amb_num(col_ndx,lyr_ndx) = amb_num(col_ndx,lyr_ndx) - delnum_sum
394 : end if
395 :
396 : ! apply the total number change to bin number
397 180451908 : cld_num(col_ndx,lyr_ndx) = cld_num(col_ndx,lyr_ndx) + delnum_sum
398 :
399 180451908 : end subroutine update_bin
400 :
401 : !------------------------------------------------------------------------------
402 : ! returns the volume-weighted fractions of aerosol subset `bin_ndx` that can act
403 : ! as heterogeneous freezing nuclei
404 : !------------------------------------------------------------------------------
405 7058880 : function hetfrz_size_wght(self, bin_ndx, ncol, nlev) result(wght)
406 : class(modal_aerosol_state), intent(in) :: self
407 : integer, intent(in) :: bin_ndx ! bin number
408 : integer, intent(in) :: ncol ! number of columns
409 : integer, intent(in) :: nlev ! number of vertical levels
410 :
411 : real(r8) :: wght(ncol,nlev)
412 :
413 : character(len=aero_name_len) :: modetype
414 :
415 5484339360 : wght(:,:) = 1._r8
416 :
417 3529440 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
418 :
419 3529440 : if (trim(modetype) == 'aitken') then
420 0 : wght(:,:) = 0._r8
421 : end if
422 :
423 3529440 : end function hetfrz_size_wght
424 :
425 : !------------------------------------------------------------------------------
426 : ! returns hygroscopicity for a given radiation diagnostic list number and
427 : ! bin number
428 : !------------------------------------------------------------------------------
429 0 : function hygroscopicity(self, list_ndx, bin_ndx) result(kappa)
430 : class(modal_aerosol_state), intent(in) :: self
431 : integer, intent(in) :: list_ndx ! rad climate list number
432 : integer, intent(in) :: bin_ndx ! bin number
433 :
434 : real(r8), pointer :: kappa(:,:) ! hygroscopicity (ncol,nlev)
435 :
436 0 : nullify(kappa)
437 :
438 0 : end function hygroscopicity
439 :
440 : !------------------------------------------------------------------------------
441 : ! returns aerosol wet diameter and aerosol water concentration for a given
442 : ! radiation diagnostic list number and bin number
443 : !------------------------------------------------------------------------------
444 990720 : subroutine water_uptake(self, aero_props, list_idx, bin_idx, ncol, nlev, dgnumwet, qaerwat)
445 : use modal_aero_wateruptake, only: modal_aero_wateruptake_dr
446 : use modal_aero_calcsize, only: modal_aero_calcsize_diag
447 :
448 : class(modal_aerosol_state), intent(in) :: self
449 : class(aerosol_properties), intent(in) :: aero_props
450 : integer, intent(in) :: list_idx ! rad climate/diags list number
451 : integer, intent(in) :: bin_idx ! bin number
452 : integer, intent(in) :: ncol ! number of columns
453 : integer, intent(in) :: nlev ! number of levels
454 : real(r8),intent(out) :: dgnumwet(ncol,nlev) ! aerosol wet diameter (m)
455 : real(r8),intent(out) :: qaerwat(ncol,nlev) ! aerosol water concentration (g/g)
456 :
457 : integer :: istat, nmodes
458 495360 : real(r8), pointer :: dgnumdry_m(:,:,:) ! number mode dry diameter for all modes
459 495360 : real(r8), pointer :: dgnumwet_m(:,:,:) ! number mode wet diameter for all modes
460 495360 : real(r8), pointer :: qaerwat_m(:,:,:) ! aerosol water (g/g) for all modes
461 495360 : real(r8), pointer :: wetdens_m(:,:,:) !
462 495360 : real(r8), pointer :: hygro_m(:,:,:) !
463 495360 : real(r8), pointer :: dryvol_m(:,:,:) !
464 495360 : real(r8), pointer :: dryrad_m(:,:,:) !
465 495360 : real(r8), pointer :: drymass_m(:,:,:) !
466 495360 : real(r8), pointer :: so4dryvol_m(:,:,:) !
467 495360 : real(r8), pointer :: naer_m(:,:,:) !
468 :
469 990720 : nmodes = aero_props%nbins()
470 :
471 495360 : if (list_idx == 0) then
472 : ! water uptake and wet radius for the climate list has already been calculated
473 495360 : call pbuf_get_field(self%pbuf, pbuf_get_index('DGNUMWET'), dgnumwet_m)
474 495360 : call pbuf_get_field(self%pbuf, pbuf_get_index('QAERWAT'), qaerwat_m)
475 :
476 769731840 : dgnumwet(:ncol,:nlev) = dgnumwet_m(:ncol,:nlev,bin_idx)
477 769731840 : qaerwat (:ncol,:nlev) = qaerwat_m(:ncol,:nlev,bin_idx)
478 :
479 : else
480 : ! If doing a diagnostic calculation then need to calculate the wet radius
481 : ! and water uptake for the diagnostic modes
482 : allocate(dgnumdry_m(ncol,nlev,nmodes), dgnumwet_m(ncol,nlev,nmodes), &
483 : qaerwat_m(ncol,nlev,nmodes), wetdens_m(ncol,nlev,nmodes), &
484 : hygro_m(ncol,nlev,nmodes), dryvol_m(ncol,nlev,nmodes), &
485 : dryrad_m(ncol,nlev,nmodes), drymass_m(ncol,nlev,nmodes), &
486 0 : so4dryvol_m(ncol,nlev,nmodes), naer_m(ncol,nlev,nmodes), stat=istat)
487 0 : if (istat > 0) then
488 0 : dgnumwet = -huge(1._r8)
489 0 : qaerwat = -huge(1._r8)
490 : return
491 : end if
492 : call modal_aero_calcsize_diag(self%state, self%pbuf, list_idx, dgnumdry_m, hygro_m, &
493 0 : dryvol_m, dryrad_m, drymass_m, so4dryvol_m, naer_m)
494 : call modal_aero_wateruptake_dr(self%state, self%pbuf, list_idx, dgnumdry_m, dgnumwet_m, &
495 : qaerwat_m, wetdens_m, hygro_m, dryvol_m, dryrad_m, &
496 0 : drymass_m, so4dryvol_m, naer_m)
497 :
498 0 : dgnumwet(:ncol,:nlev) = dgnumwet_m(:ncol,:nlev,bin_idx)
499 0 : qaerwat (:ncol,:nlev) = qaerwat_m(:ncol,:nlev,bin_idx)
500 :
501 0 : deallocate(dgnumdry_m)
502 0 : deallocate(dgnumwet_m)
503 0 : deallocate(qaerwat_m)
504 0 : deallocate(wetdens_m)
505 0 : deallocate(hygro_m)
506 0 : deallocate(dryvol_m)
507 0 : deallocate(dryrad_m)
508 0 : deallocate(drymass_m)
509 0 : deallocate(so4dryvol_m)
510 0 : deallocate(naer_m)
511 : endif
512 :
513 :
514 990720 : end subroutine water_uptake
515 :
516 : !------------------------------------------------------------------------------
517 : ! aerosol dry volume (m3/kg) for given radiation diagnostic list number and bin number
518 : !------------------------------------------------------------------------------
519 495360 : function dry_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
520 :
521 : class(modal_aerosol_state), intent(in) :: self
522 : class(aerosol_properties), intent(in) :: aero_props
523 :
524 : integer, intent(in) :: list_idx ! rad climate/diags list number
525 : integer, intent(in) :: bin_idx ! bin number
526 : integer, intent(in) :: ncol ! number of columns
527 : integer, intent(in) :: nlev ! number of levels
528 :
529 : real(r8) :: vol(ncol,nlev) ! m3/kg
530 :
531 123840 : real(r8), pointer :: mmr(:,:)
532 : real(r8) :: specdens ! species density (kg/m3)
533 :
534 : integer :: ispec
535 :
536 192432960 : vol(:,:) = 0._r8
537 :
538 588240 : do ispec = 1, aero_props%nspecies(list_idx,bin_idx)
539 464400 : call self%get_ambient_mmr(list_idx, ispec, bin_idx, mmr)
540 464400 : call aero_props%get(bin_idx, ispec, list_ndx=list_idx, density=specdens)
541 721747440 : vol(:ncol,:) = vol(:ncol,:) + mmr(:ncol,:)/specdens
542 : end do
543 :
544 619200 : end function dry_volume
545 :
546 : !------------------------------------------------------------------------------
547 : ! aerosol wet volume (m3/kg) for given radiation diagnostic list number and bin number
548 : !------------------------------------------------------------------------------
549 495360 : function wet_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
550 :
551 : class(modal_aerosol_state), intent(in) :: self
552 : class(aerosol_properties), intent(in) :: aero_props
553 :
554 : integer, intent(in) :: list_idx ! rad climate/diags list number
555 : integer, intent(in) :: bin_idx ! bin number
556 : integer, intent(in) :: ncol ! number of columns
557 : integer, intent(in) :: nlev ! number of levels
558 :
559 : real(r8) :: vol(ncol,nlev) ! m3/kg
560 :
561 123840 : real(r8) :: dryvol(ncol,nlev)
562 247680 : real(r8) :: watervol(ncol,nlev)
563 :
564 192432960 : dryvol = self%dry_volume(aero_props, list_idx, bin_idx, ncol, nlev)
565 192432960 : watervol = self%water_volume(aero_props, list_idx, bin_idx, ncol, nlev)
566 :
567 192432960 : vol = watervol + dryvol
568 :
569 123840 : end function wet_volume
570 :
571 : !------------------------------------------------------------------------------
572 : ! aerosol water volume (m3/kg) for given radiation diagnostic list number and bin number
573 : !------------------------------------------------------------------------------
574 990720 : function water_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
575 :
576 : class(modal_aerosol_state), intent(in) :: self
577 : class(aerosol_properties), intent(in) :: aero_props
578 :
579 : integer, intent(in) :: list_idx ! rad climate/diags list number
580 : integer, intent(in) :: bin_idx ! bin number
581 : integer, intent(in) :: ncol ! number of columns
582 : integer, intent(in) :: nlev ! number of levels
583 :
584 : real(r8) :: vol(ncol,nlev) ! m3/kg
585 :
586 495360 : real(r8) :: dgnumwet(ncol,nlev)
587 247680 : real(r8) :: qaerwat(ncol,nlev)
588 :
589 247680 : call self%water_uptake(aero_props, list_idx, bin_idx, ncol, nlev, dgnumwet, qaerwat)
590 :
591 384865920 : vol(:ncol,:nlev) = qaerwat(:ncol,:nlev)*rh2odens
592 384865920 : where (vol<0._r8)
593 : vol = 0._r8
594 : end where
595 :
596 247680 : end function water_volume
597 :
598 244584 : end module modal_aerosol_state_mod
|