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 : procedure :: wet_diameter
44 : procedure :: convcld_actfrac
45 :
46 : final :: destructor
47 :
48 : end type modal_aerosol_state
49 :
50 : interface modal_aerosol_state
51 : procedure :: constructor
52 : end interface modal_aerosol_state
53 :
54 : real(r8), parameter :: rh2odens = 1._r8/rhoh2o
55 :
56 : contains
57 :
58 : !------------------------------------------------------------------------------
59 : !------------------------------------------------------------------------------
60 0 : function constructor(state,pbuf) result(newobj)
61 : type(physics_state), target :: state
62 : type(physics_buffer_desc), pointer :: pbuf(:)
63 :
64 : type(modal_aerosol_state), pointer :: newobj
65 :
66 : integer :: ierr
67 :
68 0 : allocate(newobj,stat=ierr)
69 0 : if( ierr /= 0 ) then
70 0 : nullify(newobj)
71 : return
72 : end if
73 :
74 0 : newobj%state => state
75 0 : newobj%pbuf => pbuf
76 :
77 0 : end function constructor
78 :
79 : !------------------------------------------------------------------------------
80 : !------------------------------------------------------------------------------
81 0 : subroutine destructor(self)
82 : type(modal_aerosol_state), intent(inout) :: self
83 :
84 0 : nullify(self%state)
85 0 : nullify(self%pbuf)
86 :
87 0 : end subroutine destructor
88 :
89 : !------------------------------------------------------------------------------
90 : ! sets transported components
91 : ! This aerosol model with the state of the transported aerosol constituents
92 : ! (mass mixing ratios or number mixing ratios)
93 : !------------------------------------------------------------------------------
94 0 : subroutine set_transported( self, transported_array )
95 : class(modal_aerosol_state), intent(inout) :: self
96 : real(r8), intent(in) :: transported_array(:,:,:)
97 : ! to be implemented later
98 0 : end subroutine set_transported
99 :
100 : !------------------------------------------------------------------------------
101 : ! returns transported components
102 : ! This returns to current state of the transported aerosol constituents
103 : ! (mass mixing ratios or number mixing ratios)
104 : !------------------------------------------------------------------------------
105 0 : subroutine get_transported( self, transported_array )
106 : class(modal_aerosol_state), intent(in) :: self
107 : real(r8), intent(out) :: transported_array(:,:,:)
108 : ! to be implemented later
109 0 : end subroutine get_transported
110 :
111 : !------------------------------------------------------------------------
112 : ! Total aerosol mass mixing ratio for a bin in a given grid box location (column and layer)
113 : !------------------------------------------------------------------------
114 0 : function ambient_total_bin_mmr(self, aero_props, bin_ndx, col_ndx, lyr_ndx) result(mmr_tot)
115 : class(modal_aerosol_state), intent(in) :: self
116 : class(aerosol_properties), intent(in) :: aero_props
117 : integer, intent(in) :: bin_ndx ! bin index
118 : integer, intent(in) :: col_ndx ! column index
119 : integer, intent(in) :: lyr_ndx ! vertical layer index
120 :
121 : real(r8) :: mmr_tot ! mass mixing ratios totaled for all species
122 0 : real(r8),pointer :: mmrptr(:,:)
123 : integer :: spec_ndx
124 :
125 0 : mmr_tot = 0._r8
126 :
127 0 : do spec_ndx=1,aero_props%nspecies(bin_ndx)
128 0 : call rad_cnst_get_aer_mmr(0, bin_ndx, spec_ndx, 'a', self%state, self%pbuf, mmrptr)
129 0 : mmr_tot = mmr_tot + mmrptr(col_ndx,lyr_ndx)
130 : end do
131 :
132 0 : end function ambient_total_bin_mmr
133 :
134 : !------------------------------------------------------------------------------
135 : ! returns ambient aerosol mass mixing ratio for a given species index and bin index
136 : !------------------------------------------------------------------------------
137 0 : subroutine get_ambient_mmr_0list(self, species_ndx, bin_ndx, mmr)
138 : class(modal_aerosol_state), intent(in) :: self
139 : integer, intent(in) :: species_ndx ! species index
140 : integer, intent(in) :: bin_ndx ! bin index
141 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
142 :
143 0 : call rad_cnst_get_aer_mmr(0, bin_ndx, species_ndx, 'a', self%state, self%pbuf, mmr)
144 0 : end subroutine get_ambient_mmr_0list
145 :
146 : !------------------------------------------------------------------------------
147 : ! returns ambient aerosol mass mixing ratio for a given radiation diagnostics
148 : ! list index, species index and bin index
149 : !------------------------------------------------------------------------------
150 0 : subroutine get_ambient_mmr_rlist(self, list_ndx, species_ndx, bin_ndx, mmr)
151 : class(modal_aerosol_state), intent(in) :: self
152 : integer, intent(in) :: list_ndx ! rad climate list index
153 : integer, intent(in) :: species_ndx ! species index
154 : integer, intent(in) :: bin_ndx ! bin index
155 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
156 :
157 0 : call rad_cnst_get_aer_mmr(list_ndx, bin_ndx, species_ndx, 'a', self%state, self%pbuf, mmr)
158 0 : end subroutine get_ambient_mmr_rlist
159 :
160 : !------------------------------------------------------------------------------
161 : ! returns cloud-borne aerosol number mixing ratio for a given species index and bin index
162 : !------------------------------------------------------------------------------
163 0 : subroutine get_cldbrne_mmr(self, species_ndx, bin_ndx, mmr)
164 : class(modal_aerosol_state), intent(in) :: self
165 : integer, intent(in) :: species_ndx ! species index
166 : integer, intent(in) :: bin_ndx ! bin index
167 : real(r8), pointer :: mmr(:,:) ! mass mixing ratios (ncol,nlev)
168 :
169 0 : call rad_cnst_get_aer_mmr(0, bin_ndx, species_ndx, 'c', self%state, self%pbuf, mmr)
170 0 : end subroutine get_cldbrne_mmr
171 :
172 : !------------------------------------------------------------------------------
173 : ! returns ambient aerosol number mixing ratio for a given species index and bin index
174 : !------------------------------------------------------------------------------
175 0 : subroutine get_ambient_num(self, bin_ndx, num)
176 : class(modal_aerosol_state), intent(in) :: self
177 : integer, intent(in) :: bin_ndx ! bin index
178 : real(r8), pointer :: num(:,:) ! number densities
179 :
180 0 : call rad_cnst_get_mode_num(0, bin_ndx, 'a', self%state, self%pbuf, num)
181 0 : end subroutine get_ambient_num
182 :
183 : !------------------------------------------------------------------------------
184 : ! returns cloud-borne aerosol number mixing ratio for a given species index and bin index
185 : !------------------------------------------------------------------------------
186 0 : subroutine get_cldbrne_num(self, bin_ndx, num)
187 : class(modal_aerosol_state), intent(in) :: self
188 : integer, intent(in) :: bin_ndx ! bin index
189 : real(r8), pointer :: num(:,:)
190 :
191 0 : call rad_cnst_get_mode_num(0, bin_ndx, 'c', self%state, self%pbuf, num)
192 0 : end subroutine get_cldbrne_num
193 :
194 : !------------------------------------------------------------------------------
195 : ! returns interstitial and cloud-borne aerosol states
196 : !------------------------------------------------------------------------------
197 0 : subroutine get_states( self, aero_props, raer, qqcw )
198 : class(modal_aerosol_state), intent(in) :: self
199 : class(aerosol_properties), intent(in) :: aero_props
200 : type(ptr2d_t), intent(out) :: raer(:)
201 : type(ptr2d_t), intent(out) :: qqcw(:)
202 :
203 : integer :: ibin,ispc, indx
204 :
205 0 : do ibin = 1, aero_props%nbins()
206 0 : indx = aero_props%indexer(ibin, 0)
207 0 : call self%get_ambient_num(ibin, raer(indx)%fld)
208 0 : call self%get_cldbrne_num(ibin, qqcw(indx)%fld)
209 0 : do ispc = 1, aero_props%nspecies(ibin)
210 0 : indx = aero_props%indexer(ibin, ispc)
211 0 : call self%get_ambient_mmr(ispc,ibin, raer(indx)%fld)
212 0 : call self%get_cldbrne_mmr(ispc,ibin, qqcw(indx)%fld)
213 : end do
214 : end do
215 :
216 0 : end subroutine get_states
217 :
218 : !------------------------------------------------------------------------------
219 : ! return aerosol bin size weights for a given bin
220 : !------------------------------------------------------------------------------
221 0 : subroutine icenuc_size_wght_arr(self, bin_ndx, ncol, nlev, species_type, use_preexisting_ice, wght)
222 : class(modal_aerosol_state), intent(in) :: self
223 : integer, intent(in) :: bin_ndx ! bin number
224 : integer, intent(in) :: ncol ! number of columns
225 : integer, intent(in) :: nlev ! number of vertical levels
226 : character(len=*), intent(in) :: species_type ! species type
227 : logical, intent(in) :: use_preexisting_ice ! pre-existing ice flag
228 : real(r8), intent(out) :: wght(:,:)
229 :
230 : character(len=aero_name_len) :: modetype
231 0 : real(r8), pointer :: dgnum(:,:,:) ! mode dry radius
232 : real(r8) :: sigmag_aitken
233 : integer :: i,k
234 :
235 0 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
236 :
237 0 : wght = 0._r8
238 :
239 0 : select case ( trim(species_type) )
240 : case('dust')
241 0 : if (modetype=='coarse' .or. modetype=='coarse_dust') then
242 0 : wght(:ncol,:) = 1._r8
243 : end if
244 : case('sulfate')
245 0 : if (modetype=='aitken') then
246 0 : if ( use_preexisting_ice ) then
247 0 : wght(:ncol,:) = 1._r8
248 : else
249 0 : call rad_cnst_get_mode_props(0, bin_ndx, sigmag=sigmag_aitken)
250 0 : call pbuf_get_field(self%pbuf, pbuf_get_index('DGNUM' ), dgnum)
251 0 : do k = 1,nlev
252 0 : do i = 1,ncol
253 0 : if (dgnum(i,k,bin_ndx) > 0._r8) then
254 : ! only allow so4 with D>0.1 um in ice nucleation
255 0 : wght(i,k) = max(0._r8,(0.5_r8 - 0.5_r8* &
256 : erf(log(0.1e-6_r8/dgnum(i,k,bin_ndx))/ &
257 0 : (2._r8**0.5_r8*log(sigmag_aitken))) ))
258 : end if
259 : end do
260 : end do
261 : endif
262 : endif
263 : case('black-c')
264 0 : if (modetype=='accum') then
265 0 : wght(:ncol,:) = 1._r8
266 : endif
267 : case('sulfate_strat')
268 0 : if (modetype=='accum' .or. modetype=='coarse' .or. modetype=='coarse_strat') then
269 0 : wght(:ncol,:) = 1._r8
270 : endif
271 : end select
272 :
273 0 : end subroutine icenuc_size_wght_arr
274 :
275 : !------------------------------------------------------------------------------
276 : ! return aerosol bin size weights for a given bin, column and vertical layer
277 : !------------------------------------------------------------------------------
278 0 : subroutine icenuc_size_wght_val(self, bin_ndx, col_ndx, lyr_ndx, species_type, use_preexisting_ice, wght)
279 : class(modal_aerosol_state), intent(in) :: self
280 : integer, intent(in) :: bin_ndx ! bin number
281 : integer, intent(in) :: col_ndx ! column index
282 : integer, intent(in) :: lyr_ndx ! vertical layer index
283 : character(len=*), intent(in) :: species_type ! species type
284 : logical, intent(in) :: use_preexisting_ice ! pre-existing ice flag
285 : real(r8), intent(out) :: wght
286 :
287 : character(len=aero_name_len) :: modetype
288 0 : real(r8), pointer :: dgnum(:,:,:) ! mode dry radius
289 : real(r8) :: sigmag_aitken
290 :
291 0 : wght = 0._r8
292 :
293 0 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
294 :
295 0 : select case ( trim(species_type) )
296 : case('dust')
297 0 : if (modetype=='coarse' .or. modetype=='coarse_dust') then
298 0 : wght = 1._r8
299 : end if
300 : case('sulfate')
301 0 : if (modetype=='aitken') then
302 0 : if ( use_preexisting_ice ) then
303 0 : wght = 1._r8
304 : else
305 0 : call rad_cnst_get_mode_props(0, bin_ndx, sigmag=sigmag_aitken)
306 0 : call pbuf_get_field(self%pbuf, pbuf_get_index('DGNUM' ), dgnum)
307 :
308 0 : if (dgnum(col_ndx,lyr_ndx,bin_ndx) > 0._r8) then
309 : ! only allow so4 with D>0.1 um in ice nucleation
310 : wght = max(0._r8,(0.5_r8 - 0.5_r8* &
311 : erf(log(0.1e-6_r8/dgnum(col_ndx,lyr_ndx,bin_ndx))/ &
312 0 : (2._r8**0.5_r8*log(sigmag_aitken))) ))
313 :
314 : end if
315 : endif
316 : endif
317 : case('black-c')
318 0 : if (modetype=='accum') then
319 0 : wght = 1._r8
320 : endif
321 : case('sulfate_strat')
322 0 : if (modetype=='accum' .or. modetype=='coarse' .or. modetype=='coarse_strat') then
323 0 : wght = 1._r8
324 : endif
325 : end select
326 :
327 0 : end subroutine icenuc_size_wght_val
328 :
329 : !------------------------------------------------------------------------------
330 : ! returns aerosol type weights for a given aerosol type and bin
331 : !------------------------------------------------------------------------------
332 0 : subroutine icenuc_type_wght(self, bin_ndx, ncol, nlev, species_type, aero_props, rho, wght, cloud_borne)
333 :
334 : use aerosol_properties_mod, only: aerosol_properties
335 :
336 : class(modal_aerosol_state), intent(in) :: self
337 : integer, intent(in) :: bin_ndx ! bin number
338 : integer, intent(in) :: ncol ! number of columns
339 : integer, intent(in) :: nlev ! number of vertical levels
340 : character(len=*), intent(in) :: species_type ! species type
341 : class(aerosol_properties), intent(in) :: aero_props ! aerosol properties object
342 : real(r8), intent(in) :: rho(:,:) ! air density (kg m-3)
343 : real(r8), intent(out) :: wght(:,:) ! type weights
344 : logical, optional, intent(in) :: cloud_borne ! if TRUE cloud-borne aerosols are used
345 : ! otherwise ambient aerosols are used
346 :
347 : character(len=aero_name_len) :: modetype
348 :
349 0 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
350 :
351 0 : wght = 0._r8
352 :
353 0 : if (species_type == 'dust') then
354 0 : if (modetype=='coarse_dust') then
355 0 : wght(:ncol,:) = 1._r8
356 : else
357 0 : call self%icenuc_type_wght_base(bin_ndx, ncol, nlev, species_type, aero_props, rho, wght, cloud_borne)
358 : end if
359 0 : else if (species_type == 'sulfate_strat') then
360 0 : if (modetype=='accum') then
361 0 : wght(:ncol,:) = 1._r8
362 0 : elseif ( modetype=='coarse' .or. modetype=='coarse_strat') then
363 0 : call self%icenuc_type_wght_base(bin_ndx, ncol, nlev, species_type, aero_props, rho, wght, cloud_borne)
364 : endif
365 : else
366 0 : wght(:ncol,:) = 1._r8
367 : end if
368 :
369 0 : end subroutine icenuc_type_wght
370 :
371 : !------------------------------------------------------------------------------
372 : !------------------------------------------------------------------------------
373 0 : subroutine update_bin( self, bin_ndx, col_ndx, lyr_ndx, delmmr_sum, delnum_sum, tnd_ndx, dtime, tend )
374 : class(modal_aerosol_state), intent(in) :: self
375 : integer, intent(in) :: bin_ndx ! bin number
376 : integer, intent(in) :: col_ndx ! column index
377 : integer, intent(in) :: lyr_ndx ! vertical layer index
378 : real(r8),intent(in) :: delmmr_sum ! mass mixing ratio change summed over all species in bin
379 : real(r8),intent(in) :: delnum_sum ! number mixing ratio change summed over all species in bin
380 : integer, intent(in) :: tnd_ndx ! tendency index
381 : real(r8),intent(in) :: dtime ! time step size (sec)
382 : real(r8),intent(inout) :: tend(:,:,:) ! tendency
383 :
384 0 : real(r8), pointer :: amb_num(:,:)
385 0 : real(r8), pointer :: cld_num(:,:)
386 :
387 0 : call self%get_ambient_num(bin_ndx, amb_num)
388 0 : call self%get_cldbrne_num(bin_ndx, cld_num)
389 :
390 : ! if there is no bin mass compute updates/tendencies for bin number
391 : ! -- apply the total number change to bin number
392 0 : if (tnd_ndx>0) then
393 0 : tend(col_ndx,lyr_ndx,tnd_ndx) = -delnum_sum/dtime
394 : else
395 0 : amb_num(col_ndx,lyr_ndx) = amb_num(col_ndx,lyr_ndx) - delnum_sum
396 : end if
397 :
398 : ! apply the total number change to bin number
399 0 : cld_num(col_ndx,lyr_ndx) = cld_num(col_ndx,lyr_ndx) + delnum_sum
400 :
401 0 : end subroutine update_bin
402 :
403 : !------------------------------------------------------------------------------
404 : ! returns the volume-weighted fractions of aerosol subset `bin_ndx` that can act
405 : ! as heterogeneous freezing nuclei
406 : !------------------------------------------------------------------------------
407 0 : function hetfrz_size_wght(self, bin_ndx, ncol, nlev) result(wght)
408 : class(modal_aerosol_state), intent(in) :: self
409 : integer, intent(in) :: bin_ndx ! bin number
410 : integer, intent(in) :: ncol ! number of columns
411 : integer, intent(in) :: nlev ! number of vertical levels
412 :
413 : real(r8) :: wght(ncol,nlev)
414 :
415 : character(len=aero_name_len) :: modetype
416 :
417 0 : wght(:,:) = 1._r8
418 :
419 0 : call rad_cnst_get_info(0, bin_ndx, mode_type=modetype)
420 :
421 0 : if (trim(modetype) == 'aitken') then
422 0 : wght(:,:) = 0._r8
423 : end if
424 :
425 0 : end function hetfrz_size_wght
426 :
427 : !------------------------------------------------------------------------------
428 : ! returns hygroscopicity for a given radiation diagnostic list number and
429 : ! bin number
430 : !------------------------------------------------------------------------------
431 0 : function hygroscopicity(self, list_ndx, bin_ndx) result(kappa)
432 : class(modal_aerosol_state), intent(in) :: self
433 : integer, intent(in) :: list_ndx ! rad climate list number
434 : integer, intent(in) :: bin_ndx ! bin number
435 :
436 : real(r8), pointer :: kappa(:,:) ! hygroscopicity (ncol,nlev)
437 :
438 0 : nullify(kappa)
439 :
440 0 : end function hygroscopicity
441 :
442 : !------------------------------------------------------------------------------
443 : ! returns aerosol wet diameter and aerosol water concentration for a given
444 : ! radiation diagnostic list number and bin number
445 : !------------------------------------------------------------------------------
446 0 : subroutine water_uptake(self, aero_props, list_idx, bin_idx, ncol, nlev, dgnumwet, qaerwat)
447 : use modal_aero_wateruptake, only: modal_aero_wateruptake_dr
448 : use modal_aero_calcsize, only: modal_aero_calcsize_diag
449 :
450 : class(modal_aerosol_state), intent(in) :: self
451 : class(aerosol_properties), intent(in) :: aero_props
452 : integer, intent(in) :: list_idx ! rad climate/diags list number
453 : integer, intent(in) :: bin_idx ! bin number
454 : integer, intent(in) :: ncol ! number of columns
455 : integer, intent(in) :: nlev ! number of levels
456 : real(r8),intent(out) :: dgnumwet(ncol,nlev) ! aerosol wet diameter (m)
457 : real(r8),intent(out) :: qaerwat(ncol,nlev) ! aerosol water concentration (g/g)
458 :
459 : integer :: istat, nmodes
460 0 : real(r8), pointer :: dgnumdry_m(:,:,:) ! number mode dry diameter for all modes
461 0 : real(r8), pointer :: dgnumwet_m(:,:,:) ! number mode wet diameter for all modes
462 0 : real(r8), pointer :: qaerwat_m(:,:,:) ! aerosol water (g/g) for all modes
463 0 : real(r8), pointer :: wetdens_m(:,:,:) !
464 0 : real(r8), pointer :: hygro_m(:,:,:) !
465 0 : real(r8), pointer :: dryvol_m(:,:,:) !
466 0 : real(r8), pointer :: dryrad_m(:,:,:) !
467 0 : real(r8), pointer :: drymass_m(:,:,:) !
468 0 : real(r8), pointer :: so4dryvol_m(:,:,:) !
469 0 : real(r8), pointer :: naer_m(:,:,:) !
470 :
471 0 : nmodes = aero_props%nbins()
472 :
473 0 : if (list_idx == 0) then
474 : ! water uptake and wet radius for the climate list has already been calculated
475 0 : call pbuf_get_field(self%pbuf, pbuf_get_index('DGNUMWET'), dgnumwet_m)
476 0 : call pbuf_get_field(self%pbuf, pbuf_get_index('QAERWAT'), qaerwat_m)
477 :
478 0 : dgnumwet(:ncol,:nlev) = dgnumwet_m(:ncol,:nlev,bin_idx)
479 0 : qaerwat (:ncol,:nlev) = qaerwat_m(:ncol,:nlev,bin_idx)
480 :
481 : else
482 : ! If doing a diagnostic calculation then need to calculate the wet radius
483 : ! and water uptake for the diagnostic modes
484 : allocate(dgnumdry_m(ncol,nlev,nmodes), dgnumwet_m(ncol,nlev,nmodes), &
485 : qaerwat_m(ncol,nlev,nmodes), wetdens_m(ncol,nlev,nmodes), &
486 : hygro_m(ncol,nlev,nmodes), dryvol_m(ncol,nlev,nmodes), &
487 : dryrad_m(ncol,nlev,nmodes), drymass_m(ncol,nlev,nmodes), &
488 0 : so4dryvol_m(ncol,nlev,nmodes), naer_m(ncol,nlev,nmodes), stat=istat)
489 0 : if (istat > 0) then
490 0 : dgnumwet = -huge(1._r8)
491 0 : qaerwat = -huge(1._r8)
492 : return
493 : end if
494 : call modal_aero_calcsize_diag(self%state, self%pbuf, list_idx, dgnumdry_m, hygro_m, &
495 0 : dryvol_m, dryrad_m, drymass_m, so4dryvol_m, naer_m)
496 : call modal_aero_wateruptake_dr(self%state, self%pbuf, list_idx, dgnumdry_m, dgnumwet_m, &
497 : qaerwat_m, wetdens_m, hygro_m, dryvol_m, dryrad_m, &
498 0 : drymass_m, so4dryvol_m, naer_m)
499 :
500 0 : dgnumwet(:ncol,:nlev) = dgnumwet_m(:ncol,:nlev,bin_idx)
501 0 : qaerwat (:ncol,:nlev) = qaerwat_m(:ncol,:nlev,bin_idx)
502 :
503 0 : deallocate(dgnumdry_m)
504 0 : deallocate(dgnumwet_m)
505 0 : deallocate(qaerwat_m)
506 0 : deallocate(wetdens_m)
507 0 : deallocate(hygro_m)
508 0 : deallocate(dryvol_m)
509 0 : deallocate(dryrad_m)
510 0 : deallocate(drymass_m)
511 0 : deallocate(so4dryvol_m)
512 0 : deallocate(naer_m)
513 : endif
514 :
515 :
516 0 : end subroutine water_uptake
517 :
518 : !------------------------------------------------------------------------------
519 : ! aerosol dry volume (m3/kg) for given radiation diagnostic list number and bin number
520 : !------------------------------------------------------------------------------
521 0 : function dry_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
522 :
523 : class(modal_aerosol_state), intent(in) :: self
524 : class(aerosol_properties), intent(in) :: aero_props
525 :
526 : integer, intent(in) :: list_idx ! rad climate/diags list number
527 : integer, intent(in) :: bin_idx ! bin number
528 : integer, intent(in) :: ncol ! number of columns
529 : integer, intent(in) :: nlev ! number of levels
530 :
531 : real(r8) :: vol(ncol,nlev) ! m3/kg
532 :
533 0 : real(r8), pointer :: mmr(:,:)
534 : real(r8) :: specdens ! species density (kg/m3)
535 :
536 : integer :: ispec
537 :
538 0 : vol(:,:) = 0._r8
539 :
540 0 : do ispec = 1, aero_props%nspecies(list_idx,bin_idx)
541 0 : call self%get_ambient_mmr(list_idx, ispec, bin_idx, mmr)
542 0 : call aero_props%get(bin_idx, ispec, list_ndx=list_idx, density=specdens)
543 0 : vol(:ncol,:) = vol(:ncol,:) + mmr(:ncol,:)/specdens
544 : end do
545 :
546 0 : end function dry_volume
547 :
548 : !------------------------------------------------------------------------------
549 : ! aerosol wet volume (m3/kg) for given radiation diagnostic list number and bin number
550 : !------------------------------------------------------------------------------
551 0 : function wet_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
552 :
553 : class(modal_aerosol_state), intent(in) :: self
554 : class(aerosol_properties), intent(in) :: aero_props
555 :
556 : integer, intent(in) :: list_idx ! rad climate/diags list number
557 : integer, intent(in) :: bin_idx ! bin number
558 : integer, intent(in) :: ncol ! number of columns
559 : integer, intent(in) :: nlev ! number of levels
560 :
561 : real(r8) :: vol(ncol,nlev) ! m3/kg
562 :
563 0 : real(r8) :: dryvol(ncol,nlev)
564 0 : real(r8) :: watervol(ncol,nlev)
565 :
566 0 : dryvol = self%dry_volume(aero_props, list_idx, bin_idx, ncol, nlev)
567 0 : watervol = self%water_volume(aero_props, list_idx, bin_idx, ncol, nlev)
568 :
569 0 : vol = watervol + dryvol
570 :
571 0 : end function wet_volume
572 :
573 : !------------------------------------------------------------------------------
574 : ! aerosol water volume (m3/kg) for given radiation diagnostic list number and bin number
575 : !------------------------------------------------------------------------------
576 0 : function water_volume(self, aero_props, list_idx, bin_idx, ncol, nlev) result(vol)
577 :
578 : class(modal_aerosol_state), intent(in) :: self
579 : class(aerosol_properties), intent(in) :: aero_props
580 :
581 : integer, intent(in) :: list_idx ! rad climate/diags list number
582 : integer, intent(in) :: bin_idx ! bin number
583 : integer, intent(in) :: ncol ! number of columns
584 : integer, intent(in) :: nlev ! number of levels
585 :
586 : real(r8) :: vol(ncol,nlev) ! m3/kg
587 :
588 0 : real(r8) :: dgnumwet(ncol,nlev)
589 0 : real(r8) :: qaerwat(ncol,nlev)
590 :
591 0 : call self%water_uptake(aero_props, list_idx, bin_idx, ncol, nlev, dgnumwet, qaerwat)
592 :
593 0 : vol(:ncol,:nlev) = qaerwat(:ncol,:nlev)*rh2odens
594 0 : where (vol<0._r8)
595 : vol = 0._r8
596 : end where
597 :
598 0 : end function water_volume
599 :
600 : !------------------------------------------------------------------------------
601 : ! aerosol wet diameter
602 : !------------------------------------------------------------------------------
603 0 : function wet_diameter(self, bin_idx, ncol, nlev) result(diam)
604 : class(modal_aerosol_state), intent(in) :: self
605 : integer, intent(in) :: bin_idx ! bin number
606 : integer, intent(in) :: ncol ! number of columns
607 : integer, intent(in) :: nlev ! number of levels
608 :
609 : real(r8) :: diam(ncol,nlev)
610 :
611 0 : real(r8), pointer :: dgnumwet(:,:,:)
612 :
613 0 : call pbuf_get_field(self%pbuf, pbuf_get_index('DGNUMWET'), dgnumwet)
614 :
615 0 : diam(:ncol,:nlev) = dgnumwet(:ncol,:nlev,bin_idx)
616 :
617 0 : end function wet_diameter
618 :
619 : !------------------------------------------------------------------------------
620 : ! prescribed aerosol activation fraction for convective cloud
621 : !------------------------------------------------------------------------------
622 0 : function convcld_actfrac(self, ibin, ispc, ncol, nlev) result(frac)
623 :
624 : use modal_aero_data
625 :
626 : class(modal_aerosol_state), intent(in) :: self
627 : integer, intent(in) :: ibin ! bin index
628 : integer, intent(in) :: ispc ! species index
629 : integer, intent(in) :: ncol ! number of columns
630 : integer, intent(in) :: nlev ! number of vertical levels
631 :
632 : real(r8) :: frac(ncol,nlev)
633 :
634 0 : real(r8) :: f_act_conv_coarse(ncol,nlev)
635 : real(r8) :: f_act_conv_coarse_dust, f_act_conv_coarse_nacl
636 : real(r8) :: tmpdust, tmpnacl
637 : integer :: lcoardust, lcoarnacl
638 : integer :: i,k
639 :
640 0 : f_act_conv_coarse(:,:) = 0.60_r8
641 0 : f_act_conv_coarse_dust = 0.40_r8
642 0 : f_act_conv_coarse_nacl = 0.80_r8
643 0 : if (modeptr_coarse > 0) then
644 0 : lcoardust = lptr_dust_a_amode(modeptr_coarse)
645 0 : lcoarnacl = lptr_nacl_a_amode(modeptr_coarse)
646 0 : if ((lcoardust > 0) .and. (lcoarnacl > 0)) then
647 0 : do k = 1, nlev
648 0 : do i = 1, ncol
649 0 : tmpdust = max( 0.0_r8, self%state%q(i,k,lcoardust) )
650 0 : tmpnacl = max( 0.0_r8, self%state%q(i,k,lcoarnacl) )
651 0 : if ((tmpdust+tmpnacl) > 1.0e-30_r8) then
652 0 : f_act_conv_coarse(i,k) = (f_act_conv_coarse_dust*tmpdust &
653 0 : + f_act_conv_coarse_nacl*tmpnacl)/(tmpdust+tmpnacl)
654 : end if
655 : end do
656 : end do
657 : end if
658 : end if
659 :
660 0 : if (ibin == modeptr_pcarbon) then
661 0 : frac = 0.0_r8
662 0 : else if ((ibin == modeptr_finedust) .or. (ibin == modeptr_coardust)) then
663 0 : frac = 0.4_r8
664 : else
665 0 : frac = 0.8_r8
666 : end if
667 :
668 : ! set f_act_conv for interstitial (lphase=1) coarse mode species
669 : ! for the convective in-cloud, we conceptually treat the coarse dust and seasalt
670 : ! as being externally mixed, and apply f_act_conv = f_act_conv_coarse_dust/nacl to dust/seasalt
671 : ! number and sulfate are conceptually partitioned to the dust and seasalt
672 : ! on a mass basis, so the f_act_conv for number and sulfate are
673 : ! mass-weighted averages of the values used for dust/seasalt
674 0 : if (ibin == modeptr_coarse) then
675 0 : frac = f_act_conv_coarse
676 0 : if (ispc>0) then
677 0 : if (lmassptr_amode(ispc,ibin) == lptr_dust_a_amode(ibin)) then
678 0 : frac = f_act_conv_coarse_dust
679 0 : else if (lmassptr_amode(ispc,ibin) == lptr_nacl_a_amode(ibin)) then
680 0 : frac = f_act_conv_coarse_nacl
681 : end if
682 : end if
683 : end if
684 :
685 0 : end function convcld_actfrac
686 :
687 0 : end module modal_aerosol_state_mod
|