CAMP 1.0.0
Chemistry Across Multiple Phases
rxn_HL_phase_transfer.F90
Go to the documentation of this file.
1! Copyright (C) 2021 Barcelona Supercomputing Center and University of
2! Illinois at Urbana-Champaign
3! SPDX-License-Identifier: MIT
4
5!> \file
6!> The camp_rxn_HL_phase_transfer module.
7
8!> \page camp_rxn_HL_phase_transfer CAMP: Henry's Law Phase-Transfer Reaction
9!!
10!! Henry's Law phase-trasfer reactions use equilibrium rate constants that
11!! are calculated as:
12!!
13!! \f[
14!! HLC(298K) * e^{C({1/T-1/298})}
15!! \f]
16!!
17!! where \f$HLC(298K)\f$ is the Henry's Law constant at 298 K
18!! [M\f$\mbox{Pa}^{-1}\f$], \f$C\f$ is a constant [K] and
19!! \f$T\f$ is the temperature [K]. Uptake kinetics are based on
20!! the particle effective radius, \f$r_{eff}\f$ [\f$\mbox{m}\f$], the
21!! condensing species gas-phase diffusion coefficient, \f$D_g\f$
22!! [\f$\mbox{m}^2\,\mbox{s}^{-1}\f$], its molecular weight \f$MW\f$
23!! [kg\f$\mbox{mol}^{-3}\f$], and \f$N^{*}\f$, which is
24!! used to calculate the mass accomodation coefficient.
25!!
26!! Mass accomodation coefficients and condensation rate constants are
27!! calculated using the method of Ervans et al. (2003) \cite Ervens2003 and
28!! references therein. Mass accomodation coefficients (\f$\alpha\f$) are
29!! calculated as:
30!!
31!! \f[
32!! \Delta H_{obs} = -10 \times (N^*-1) + 7.53 \times (N^{*2/3}-1) - 0.1 \times 10 \quad (\mbox{kcal}\,\mbox{M}^{-1})
33!! \f]
34!! \f[
35!! \Delta S_{obs} = -13 \times (N^*-1) - 19 \times (N^*-1) + 9.21 \times (N^{*2/3}-1) - 0.1 \times 13 \quad (\mbox{cal}\,\mbox{M}^{-1}\,\mbox{K}^{-1})
36!! \f]
37!! \f[
38!! \frac{\alpha}{1-\alpha} = e^{\frac{-\Delta G^{*}}{RT}}
39!! \f]
40!! If \f$\Delta H\f$ and \f$\Delta S\f$ are not provided, \f$\alpha\f$ is
41!! set to 0.1 \cite Zaveri2008.
42!!
43!! Condensation rate constants are calculated following \cite Zaveri2008 as:
44!! \f[
45!! k_c = 4 \pi r D_g f_{fs}( K_n, \alpha )
46!! \f]
47!! where \f$r\f$ is the radius of the particle(s) [m], \f$D_g\f$ is the
48!! diffusion coefficient of the gas-phase species [\f$\mbox{m}^2\,\mbox{s}^{-1}\f$],
49!! \f$f_{fs}( K_n, \alpha )\f$ is the Fuchs Sutugin transition regime correction
50!! factor [unitless], \f$K_n\f$ is the Knudsen Number [unitess], and \f$\alpha\f$
51!! is the mass accomodation coefficient.
52!!
53!! Rates can be calculated as:
54!! \f[
55!! r_c = [G] N_a k_c
56!! \f]
57!! where \f$[G]\f$ is the gas-phase species concentration [ppm], \f$N_a\f$ is
58!! the number concentration of particles [\f$\mbox{particle}\,\mbox{m}^{-3}\f$] and
59!! the rate \f$r_c\f$ is in [\f$\mbox{ppm}\,\mbox{s}^{-1}\f$].
60!! The particle radius used
61!! to calculate \f$k_{f}\f$ is the effective radius [\f$r_{eff}\f$], which is
62!! taken as the "least-wrong" choice for condensation rates, as it is weighted
63!! to surface area \cite Zender2002 .
64!!
65!! Input data for Phase transfer equations have the following format :
66!! \code{.json}
67!! {
68!! "type" : "HL_PHASE_TRANSFER",
69!! "gas-phase species" : "my gas spec",
70!! "aerosol-phase species" : "my aero spec",
71!! "areosol phase" : "my aqueous phase",
72!! "aerosol-phase water" : "H2O_aq",
73!! ...
74!! }
75!! \endcode
76!! The key-value pairs \b gas-phase \b species, and \b aerosol-phase
77!! \b species are required. Only one gas- and one aerosol-phase species are
78!! allowed per phase-transfer reaction. Additionally, gas-phase species must
79!! include parameters named \b HLC(298K) \b [\b M \b Pa-1], which is the Henry's
80!! Law constant at 298 K, \b HLC \b exp \b factor \b [\b K], which is the
81!! Henry's Law constant exponential factor "C", \b diffusion \b coeff \b [\b m2
82!! \b s-1], which specifies the diffusion coefficient in
83!! \f$\mbox{m}^2\,\mbox{s}^{-1}\f$, and \b molecular \b weight
84!! \b [\b kg \b mol-1], which specifies the molecular weight of the species in
85!! \f$\mbox{kg}\,\mbox{mol}^{-1}\f$. They may optionally include the
86!! parameter \b N \b star, which will be used to calculate the mass
87!! accomodation coefficient. When this parameter is not included, the mass
88!! accomodation coefficient is assumed to be 1.0.
89
90!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
91
92!> The rxn_HL_phase_transfer_t type and associated functions.
94
98 use camp_constants, only: const
100 use camp_property
101 use camp_rxn_data
102 use camp_util, only: i_kind, dp, to_string, &
105
106 implicit none
107 private
108
109#define DELTA_H_ this%condensed_data_real(1)
110#define DELTA_S_ this%condensed_data_real(2)
111#define DIFF_COEFF_ this%condensed_data_real(3)
112#define PRE_C_AVG_ this%condensed_data_real(4)
113#define A_ this%condensed_data_real(5)
114#define C_ this%condensed_data_real(6)
115#define CONV_ this%condensed_data_real(7)
116#define MW_ this%condensed_data_real(8)
117#define NUM_AERO_PHASE_ this%condensed_data_int(1)
118#define GAS_SPEC_ this%condensed_data_int(2)
119#define NUM_INT_PROP_ 2
120#define NUM_REAL_PROP_ 8
121#define NUM_ENV_PARAM_ 4
122#define DERIV_ID_(x) this%condensed_data_int(NUM_INT_PROP_+x)
123#define JAC_ID_(x) this%condensed_data_int(NUM_INT_PROP_+1+NUM_AERO_PHASE_+x)
124#define PHASE_INT_LOC_(x) this%condensed_data_int(NUM_INT_PROP_+2+6*NUM_AERO_PHASE_+x)
125#define PHASE_REAL_LOC_(x) this%condensed_data_int(NUM_INT_PROP_+2+7*NUM_AERO_PHASE_+x)
126#define AERO_SPEC_(x) this%condensed_data_int(PHASE_INT_LOC_(x))
127#define AERO_WATER_(x) this%condensed_data_int(PHASE_INT_LOC_(x)+1)
128#define AERO_PHASE_ID_(x) this%condensed_data_int(PHASE_INT_LOC_(x)+2)
129#define AERO_REP_ID_(x) this%condensed_data_int(PHASE_INT_LOC_(x)+3)
130#define NUM_AERO_PHASE_JAC_ELEM_(x) this%condensed_data_int(PHASE_INT_LOC_(x)+4)
131#define PHASE_JAC_ID_(x,s,e) this%condensed_data_int(PHASE_INT_LOC_(x)+4+(s-1)*NUM_AERO_PHASE_JAC_ELEM_(x)+e)
132#define SMALL_WATER_CONC_(x) this%condensed_data_real(PHASE_REAL_LOC_(x))
133#define EFF_RAD_JAC_ELEM_(x,e) this%condensed_data_real(PHASE_REAL_LOC_(x)+e)
134#define NUM_CONC_JAC_ELEM_(x,e) this%condensed_data_real(PHASE_REAL_LOC_(x)+NUM_AERO_PHASE_JAC_ELEM_(x)+e)
135
137
138 !> Generic test reaction data type
139 type, extends(rxn_data_t) :: rxn_hl_phase_transfer_t
140 contains
141 !> Reaction initialization
142 procedure :: initialize
143 !> Finalize the reaction
146
147 !> Constructor for rxn_HL_phase_transfer_t
149 procedure :: constructor
150 end interface rxn_hl_phase_transfer_t
151
152contains
153
154!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
155
156 !> Constructor for Phase transfer reaction
157 function constructor() result(new_obj)
158
159 !> A new reaction instance
160 type(rxn_hl_phase_transfer_t), pointer :: new_obj
161
162 allocate(new_obj)
163 new_obj%rxn_phase = aero_rxn
164
165 end function constructor
166
167!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
168
169 !> Initialize the reaction data, validating component data and loading
170 !! any required information into the condensed data arrays for use during
171 !! solving
172 subroutine initialize(this, chem_spec_data, aero_phase, aero_rep, n_cells)
173
174 !> Reaction data
175 class(rxn_hl_phase_transfer_t), intent(inout) :: this
176 !> Chemical species data
177 type(chem_spec_data_t), intent(in) :: chem_spec_data
178 !> Aerosol phase data
179 type(aero_phase_data_ptr), intent(in) :: aero_phase(:)
180 !> Aerosol representations
181 type(aero_rep_data_ptr), pointer, intent(in) :: aero_rep(:)
182 !> Number of grid cells to solve simultaneously
183 integer(kind=i_kind), intent(in) :: n_cells
184
185 type(property_t), pointer :: spec_props
186 character(len=:), allocatable :: key_name, gas_spec_name, &
187 water_name, aero_spec_name, phase_name, error_msg
188 integer(kind=i_kind) :: i_spec, i_aero_rep, n_aero_ids, i_aero_id, &
189 n_aero_jac_elem, i_phase, tmp_size
190 type(string_t), allocatable :: unique_spec_names(:), &
191 unique_water_names(:)
192 integer(kind=i_kind), allocatable :: phase_ids(:)
193 real(kind=dp) :: temp_real, n_star
194
195 ! Get the property set
196 if (.not. associated(this%property_set)) call die_msg(318525776, &
197 "Missing property set needed to initialize reaction")
198
199 ! Get the gas-phase species name
200 key_name = "gas-phase species"
201 call assert_msg(847983010, &
202 this%property_set%get_string(key_name, gas_spec_name), &
203 "Missing gas-phase species in phase-transfer reaction")
204
205 ! Get the aerosol phase name
206 key_name = "aerosol phase"
207 call assert_msg(448087197, &
208 this%property_set%get_string(key_name, phase_name), &
209 "Missing aerosol phase in phase-transfer reaction")
210
211 ! Get the aerosol-phase species name
212 key_name = "aerosol-phase species"
213 call assert_msg(797902545, &
214 this%property_set%get_string(key_name, aero_spec_name), &
215 "Missing aerosol-phase species in phase-transfer reaction")
216
217 ! Set up a general error message
218 error_msg = " for HL partitioning reaction of gas-phase species '"// &
219 gas_spec_name//"' to aerosol-phase species '"// &
220 aero_spec_name//"' in phase '"//phase_name
221
222 ! Get the aerosol-phase water name
223 key_name = "aerosol-phase water"
224 call assert_msg(386889865, &
225 this%property_set%get_string(key_name, water_name), &
226 "Missing aerosol-phase water"//error_msg)
227
228 ! Check for aerosol representations
229 call assert_msg(234155350, associated(aero_rep), &
230 "Missing aerosol representation"//error_msg)
231 call assert_msg(207961800, size(aero_rep).gt.0, &
232 "Missing aerosol representation"//error_msg)
233
234 ! Count the instances of this phase/species pair, and the number of
235 ! Jacobian elements needed in calculations of mass, volume, etc.
236 n_aero_ids = 0
237 n_aero_jac_elem = 0
238 do i_aero_rep = 1, size(aero_rep)
239
240 ! Get the unique names in this aerosol representation for the
241 ! partitioning species and aerosol-phase water
242 unique_spec_names = aero_rep(i_aero_rep)%val%unique_names( &
243 phase_name = phase_name, spec_name = aero_spec_name, &
244 phase_is_at_surface = .true.)
245 unique_water_names = aero_rep(i_aero_rep)%val%unique_names( &
246 phase_name = phase_name, spec_name = water_name, &
247 phase_is_at_surface = .true.)
248
249 ! Skip aerosol representations that do not contain this phase
250 if (.not.allocated(unique_spec_names)) cycle
251
252 ! Check the size of the unique name lists
253 call assert_msg(598091463, size(unique_spec_names).eq. &
254 size(unique_water_names), "Missing species "// &
255 aero_spec_name//" or "//water_name//" in phase "//phase_name// &
256 " or improper implementation of aerosol phase in aerosol "// &
257 "representation"//error_msg)
258
259 ! Add these instances to the list
260 n_aero_ids = n_aero_ids + size(unique_spec_names)
261
262 ! Get the number of Jacobian elements for calculations of mass, volume,
263 ! number, etc. for this partitioning into this phase
264 phase_ids = aero_rep(i_aero_rep)%val%phase_ids(phase_name, is_at_surface=.true.)
265 do i_phase = 1, size(phase_ids)
266 n_aero_jac_elem = n_aero_jac_elem + &
267 aero_rep(i_aero_rep)%val%num_jac_elem(phase_ids(i_phase))
268 end do
269
270 deallocate(unique_spec_names)
271 deallocate(unique_water_names)
272
273 end do
274
275 ! Allocate space in the condensed data arrays
276 allocate(this%condensed_data_int(num_int_prop_ + 2 + n_aero_ids * 13 + &
277 n_aero_jac_elem * 2))
278 allocate(this%condensed_data_real(num_real_prop_ + n_aero_ids + &
279 n_aero_jac_elem * 2))
280 this%condensed_data_int(:) = int(0, kind=i_kind)
281 this%condensed_data_real(:) = real(0.0, kind=dp)
282
283 ! Save space for the environment-dependent parameters
284 this%num_env_params = num_env_param_
285
286 ! Set the number of aerosol-species instances
287 num_aero_phase_ = n_aero_ids
288
289 ! Get the properties required of the aerosol species
290 call assert_msg(669162256, &
291 chem_spec_data%get_property_set(aero_spec_name, spec_props), &
292 "Missing aerosol species properties"//error_msg)
293
294 ! Get the aerosol species molecular weight
295 key_name = "molecular weight [kg mol-1]"
296 call assert_msg(209812557, spec_props%get_real(key_name, mw_), &
297 "Missing property 'MW' for aerosol species"//error_msg)
298
299 ! Set the kg/m3 -> ppm conversion prefactor (multiply by T/P to get
300 ! conversion)
301 conv_ = const%univ_gas_const / mw_ * 1.0e6
302
303 ! Get the aerosol-phase water species
304 key_name = "aerosol-phase water"
305 call assert_msg(374667967, &
306 this%property_set%get_string(key_name, water_name), &
307 "Missing aerosol-phase water"//error_msg)
308
309 ! Set the ids of each aerosol-phase species instance
310 i_aero_id = 1
311 phase_int_loc_(i_aero_id) = num_int_prop_+8*num_aero_phase_+3
312 phase_real_loc_(i_aero_id) = num_real_prop_+1
313 do i_aero_rep = 1, size(aero_rep)
314
315 ! Get the unique names in this aerosol representation for the
316 ! partitioning species and aerosol-phase water
317 unique_spec_names = aero_rep(i_aero_rep)%val%unique_names( &
318 phase_name = phase_name, spec_name = aero_spec_name, &
319 phase_is_at_surface = .true.)
320 unique_water_names = aero_rep(i_aero_rep)%val%unique_names( &
321 phase_name = phase_name, spec_name = water_name, &
322 phase_is_at_surface = .true.)
323
324 ! Get the phase ids for this aerosol phase
325 phase_ids = aero_rep(i_aero_rep)%val%phase_ids(phase_name, is_at_surface=.true.)
326
327 ! Add the species concentration and activity coefficient ids to
328 ! the condensed data, and set the number of jacobian elements for
329 ! the aerosol representations and the locations of the real data
330 do i_spec = 1, size(unique_spec_names)
331 num_aero_phase_jac_elem_(i_aero_id) = &
332 aero_rep(i_aero_rep)%val%num_jac_elem(phase_ids(i_spec))
333 aero_spec_(i_aero_id) = &
334 aero_rep(i_aero_rep)%val%spec_state_id( &
335 unique_spec_names(i_spec)%string)
336 aero_water_(i_aero_id) = &
337 aero_rep(i_aero_rep)%val%spec_state_id( &
338 unique_water_names(i_spec)%string)
339 aero_phase_id_(i_aero_id) = phase_ids(i_spec)
340 aero_rep_id_(i_aero_id) = i_aero_rep
341 i_aero_id = i_aero_id + 1
342 if (i_aero_id .le. num_aero_phase_) then
343 phase_int_loc_(i_aero_id) = phase_int_loc_(i_aero_id - 1) + 5 + &
344 2*num_aero_phase_jac_elem_(i_aero_id - 1)
345 phase_real_loc_(i_aero_id) = phase_real_loc_(i_aero_id - 1) + 1 + &
346 2*num_aero_phase_jac_elem_(i_aero_id - 1)
347 end if
348 end do
349
350 deallocate(unique_spec_names)
351 deallocate(unique_water_names)
352
353 end do
354
355 ! Save the index of the gas-phase species in the state variable array
356 gas_spec_ = chem_spec_data%gas_state_id(gas_spec_name)
357
358 ! Make sure the species exists
359 call assert_msg(955306778, gas_spec_.gt.0, &
360 "Missing gas-phase species"//error_msg)
361
362 ! Get the required properties for the gas-phase species
363 call assert_msg(757296139, &
364 chem_spec_data%get_property_set(gas_spec_name, spec_props), &
365 "Missing gas-phase species properties"//error_msg)
366
367 ! Get Henry's Law constant parameters
368 key_name = "HLC(298K) [M Pa-1]"
369 call assert_msg(637925661, spec_props%get_real(key_name, a_), &
370 "Missing Henry's Law constant at 298 K"//error_msg)
371
372 key_name = "HLC exp factor [K]"
373 call assert_msg(801365019, spec_props%get_real(key_name, c_), &
374 "Missing Henry's Law constant exponential factor"// &
375 error_msg)
376
377 ! Get N* to calculate the mass accomodation coefficient. If it is not
378 ! present, set DELTA_H_ and DELTA_S_ to zero to indicate a mass accomodation
379 ! coefficient of 1.0
380 ! Mass accomodation equation is based on equations in:
381 ! Ervens, B., et al., 2003. "CAPRAM 2.4 (MODAC mechanism): An extended
382 ! and condensed tropospheric aqueous mechanism and its application."
383 ! J. Geophys. Res. 108, 4426. doi:10.1029/2002JD002202
384 key_name = "N star"
385 if (spec_props%get_real(key_name, n_star)) then
386 ! enthalpy change (kcal mol-1)
387 delta_h_ = real(- 10.0d0*(n_star-1.0d0) + &
388 7.53d0*(n_star**(2.0d0/3.0d0)-1.0d0) - 1.0d0, kind=dp)
389 ! entropy change (cal mol-1)
390 delta_s_ = real(- 32.0d0*(n_star-1.0d0) + &
391 9.21d0*(n_star**(2.0d0/3.0d0)-1.0d0) - 1.3d0, kind=dp)
392 ! Convert dH and dS to (J mol-1)
393 delta_h_ = real(delta_h_ * 4184.0d0, kind=dp)
394 delta_s_ = real(delta_s_ * 4.184d0, kind=dp)
395 else
396 delta_h_ = real(0.0, kind=dp)
397 delta_s_ = real(0.0, kind=dp)
398 end if
399
400 ! Get the diffusion coefficient (m^2/s)
401 key_name = "diffusion coeff [m2 s-1]"
402 call assert_msg(100205531, spec_props%get_real(key_name, diff_coeff_), &
403 "Missing diffusion coefficient for gas-phase species"//error_msg)
404
405 ! Calculate the constant portion of c_rms [m/(K^2*s)]
406 key_name = "molecular weight [kg mol-1]"
407 call assert_msg(469582180, spec_props%get_real(key_name, temp_real), &
408 "Missing molecular weight for gas-phase species"//error_msg)
409 pre_c_avg_ = sqrt(8.0*const%univ_gas_const/(const%pi*temp_real))
410
411 ! Check the sizes of the data arrays
412 tmp_size = phase_int_loc_(i_aero_id - 1) + 5 + &
413 2*num_aero_phase_jac_elem_(i_aero_id - 1) - 1
414 call assert_msg(881234422, size(this%condensed_data_int) .eq. tmp_size, &
415 "int array size mismatch"//error_msg)
416 tmp_size = phase_real_loc_(i_aero_id - 1) + 1 + &
417 2*num_aero_phase_jac_elem_(i_aero_id - 1) - 1
418 call assert_msg(520767976, size(this%condensed_data_real) .eq. tmp_size,&
419 "real array size mismatch"//error_msg)
420
421 end subroutine initialize
422
423!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
424
425 !> Finalize the reaction
426 subroutine finalize(this)
427
428 !> Reaction data
429 type(rxn_hl_phase_transfer_t), intent(inout) :: this
430
431 if (associated(this%property_set)) &
432 deallocate(this%property_set)
433 if (allocated(this%condensed_data_real)) &
434 deallocate(this%condensed_data_real)
435 if (allocated(this%condensed_data_int)) &
436 deallocate(this%condensed_data_int)
437
438 end subroutine finalize
439
440!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
441
442 !> Finalize an array of reactions
443 subroutine finalize_array(this)
444
445 !> Array of reaction data
446 type(rxn_hl_phase_transfer_t), intent(inout) :: this(:)
447
448 integer(kind=i_kind) :: i
449
450 do i = 1, size(this)
451 call finalize(this(i))
452 end do
453
454 end subroutine finalize_array
455
456!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
457
Initialize the aerosol representation data, validating component data and loading any required inform...
Get the non-unique name of a chemical species by its unique name.
Interface for to_string functions.
Definition util.F90:32
The abstract aero_phase_data_t structure and associated subroutines.
subroutine finalize_array(this)
Finalize the aerosol phase data.
type(aero_phase_data_t) function, pointer constructor(phase_name, init_size)
Constructor for aero_phase_data_t.
subroutine finalize(this)
Finalize the aerosol phase data.
The abstract aero_rep_data_t structure and associated subroutines.
The camp_state_t structure and associated subroutines.
Definition camp_state.F90:9
The chem_spec_data_t structure and associated subroutines.
Physical constants.
Definition constants.F90:9
integer, parameter dp
Kind of a double precision real number.
Definition constants.F90:16
type(const_t), save const
Fixed variable for accessing the constant's values.
Definition constants.F90:77
integer, parameter i_kind
Kind of an integer.
Definition constants.F90:21
The property_t structure and associated subroutines.
Definition property.F90:9
The rxn_data_t structure and associated subroutines.
Definition rxn_data.F90:60
integer(kind=i_kind), parameter, public aero_rxn
Aerosol-phase reaction.
Definition rxn_data.F90:89
The rxn_HL_phase_transfer_t type and associated functions.
Common utility subroutines.
Definition util.F90:9
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition util.F90:165
subroutine die_msg(code, error_msg)
Error immediately.
Definition util.F90:196
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition util.F90:130
Pointer type for building arrays.
Pointer to aero_rep_data_t extending types.
Abstract reaction data type.
Definition rxn_data.F90:99
String type for building arrays of string of various size.
Definition util.F90:53