CAMP 1.0.0
Chemistry Across Multiple Phases
rxn_SIMPOL_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_SIMPOL_phase_transfer module.
7
8!> \page camp_rxn_SIMPOL_phase_transfer CAMP: SIMPOL.1 Phase-Transfer Reaction
9!!
10!! SIMPOL phase transfer reactions are based on the SIMPOL.1 model
11!! calculations of vapor pressure described by Pankow and Asher (2008)
12!! \cite Pankow2008.
13!!
14!! Mass accomodation coefficients and condensation rate constants are
15!! calculated using the method of Ervans et al. (2003) \cite Ervens2003 and
16!! references therein. Mass accomodation coefficients (\f$\alpha\f$) are
17!! calculated as:
18!!
19!! \f[
20!! \Delta H_{obs} = -10 \times (N^*-1) + 7.53 \times (N^{*2/3}-1) - 0.1 \times 10 \quad (\mbox{kcal}\,\mbox{M}^{-1})
21!! \f]
22!! \f[
23!! \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})
24!! \f]
25!! \f[
26!! \frac{\alpha}{1-\alpha} = e^{\frac{-\Delta G^{*}}{RT}}
27!! \f]
28!! If \f$\Delta H\f$ and \f$\Delta S\f$ are not provided, the mass accomodation
29!! coefficient is assumed to be 0.1 (\cite Zaveri2008).
30!!
31!! Condensation rate constants are calculated following \cite Zaveri2008 as:
32!! \f[
33!! k_c = 4 \pi r D_g f_{fs}( K_n, \alpha )
34!! \f]
35!! where \f$r\f$ is the radius of the particle(s) [m], \f$D_g\f$ is the
36!! diffusion coefficient of the gas-phase species [\f$\mbox{m}^2\,\mbox{s}^{-1}\f$],
37!! \f$f_{fs}( K_n, \alpha )\f$ is the Fuchs Sutugin transition regime correction
38!! factor [unitless], \f$K_n\f$ is the Knudsen Number [unitess], and \f$\alpha\f$
39!! is the mass accomodation coefficient.
40!!
41!! Rates can be calculated as:
42!! \f[
43!! r_c = [G] N_a k_c
44!! \f]
45!! where \f$[G]\f$ is the gas-phase species concentration [ppm], \f$N_a\f$ is
46!! the number concentration of particles [\f$\mbox{particle}\,\mbox{m}^{-3}\f$] and
47!! the rate \f$r_c\f$ is in [\f$\mbox{ppm}\,\mbox{s}^{-1}\f$].
48!! The particle radius used
49!! to calculate \f$k_{f}\f$ is the effective radius [\f$r_{eff}\f$], which is
50!! taken as the "least-wrong" choice for condensation rates, as it is weighted
51!! to surface area \cite Zender2002 .
52!!
53!! Input data for SIMPOL phase transfer reactions have the following format :
54!! \code{.json}
55!! {
56!! "type" : "SIMPOL_PHASE_TRANSFER",
57!! "gas-phase species" : "my gas spec",
58!! "aerosol phase" : "my aero phase",
59!! "aerosol-phase species" : "my aero spec",
60!! "aerosol-phase activity coefficient" : "my aero act coeff"
61!! "B" : [ 123.2e3, -41.24, 2951.2, -1.245e-4 ]
62!! ...
63!! }
64!! \endcode
65!! The key-value pairs \b gas-phase \b species, \b aerosol \b phase and
66!! \b aerosol-phase \b species are required. Only one gas- and one
67!! aerosol-phase species are allowed per phase-transfer reaction. The
68!! key-value pair \b aerosol-phase \b activity \b coefficient is optional.
69!! When it is included its value must be the name of an species of type
70!! \c ACTIVITY_COEFF that is present in the specified aerosol phase. When
71!! it is not included, activity coefficients are assume to be 1.0.
72!!
73!! Gas-phase species must include parameters named
74!! \b diffusion \b coeff \b [\b m2 \b s-1], which specifies the diffusion
75!! coefficient in \f$\mbox{m}^2\,\mbox{s}^{-1}\f$, and \b molecular
76!! \b weight \b [\b kg \b mol-1], which specifies the molecular weight of the
77!! species in \f$\mbox{kg}\,\mbox{mol}^{-1}\f$. They may optionally
78!! include the parameter \b N \b star, which will be used to calculate th
79!! mass accomodation coefficient. When this parameter is not included, the
80!! mass accomodation coefficient is assumed to be 1.0.
81!!
82!! The key-value pair \b B is also required and must have a value of an array
83!! of exactly four members that specifies the SIMPOL parameters for the
84!! partitioning species. The \b B parameters can be obtained by summing the
85!! contributions of each functional group present in the partitioning species
86!! to the overall \f$B_{n,i}\f$ for species \f$i\f$, such that:
87!! \f[
88!! B_{n,i} = \sum_{k} \nu_{k,i} B_{n,k} \forall n \in [1...4]
89!! \f]
90!! where \f$\nu_{k,i}\f$ is the number of functional groups \f$k\f$ in species
91!! \f$i\f$ and the parameters \f$B_{n,k}\f$ for each functional group \f$k\f$
92!! can be found in table 5 of Pankow and Asher (2008) \cite Pankow2008.
93
94!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
95
96!> The rxn_SIMPOL_phase_transfer_t type and associated functions.
98
102 use camp_constants, only: const
104 use camp_property
105 use camp_rxn_data
106 use camp_util, only: i_kind, dp, to_string, &
109
110 implicit none
111 private
112
113#define DELTA_H_ this%condensed_data_real(1)
114#define DELTA_S_ this%condensed_data_real(2)
115#define DIFF_COEFF_ this%condensed_data_real(3)
116#define PRE_C_AVG_ this%condensed_data_real(4)
117#define B1_ this%condensed_data_real(5)
118#define B2_ this%condensed_data_real(6)
119#define B3_ this%condensed_data_real(7)
120#define B4_ this%condensed_data_real(8)
121#define CONV_ this%condensed_data_real(9)
122#define MW_ this%condensed_data_real(10)
123#define NUM_AERO_PHASE_ this%condensed_data_int(1)
124#define GAS_SPEC_ this%condensed_data_int(2)
125#define NUM_INT_PROP_ 2
126#define NUM_REAL_PROP_ 10
127#define NUM_ENV_PARAM_ 4
128#define AERO_SPEC_(x) this%condensed_data_int(NUM_INT_PROP_+x)
129#define AERO_ACT_ID_(x) this%condensed_data_int(NUM_INT_PROP_+NUM_AERO_PHASE_+x)
130#define AERO_PHASE_ID_(x) this%condensed_data_int(NUM_INT_PROP_+2*NUM_AERO_PHASE_+x)
131#define AERO_REP_ID_(x) this%condensed_data_int(NUM_INT_PROP_+3*NUM_AERO_PHASE_+x)
132#define DERIV_ID_(x) this%condensed_data_int(NUM_INT_PROP_+4*NUM_AERO_PHASE_+x)
133#define GAS_ACT_JAC_ID_(x) this%condensed_data_int(NUM_INT_PROP_+1+5*NUM_AERO_PHASE_+x)
134#define AERO_ACT_JAC_ID_(x) this%condensed_data_int(NUM_INT_PROP_+1+6*NUM_AERO_PHASE_+x)
135#define JAC_ID_(x) this%condensed_data_int(NUM_INT_PROP_+1+7*NUM_AERO_PHASE_+x)
136#define PHASE_INT_LOC_(x) this%condensed_data_int(NUM_INT_PROP_+2+10*NUM_AERO_PHASE_+x)
137#define PHASE_REAL_LOC_(x) this%condensed_data_int(NUM_INT_PROP_+2+11*NUM_AERO_PHASE_+x)
138#define NUM_AERO_PHASE_JAC_ELEM_(x) this%condensed_data_int(PHASE_INT_LOC_(x))
139#define PHASE_JAC_ID_(x,s,e) this%condensed_data_int(PHASE_INT_LOC_(x)+(s-1)*NUM_AERO_PHASE_JAC_ELEM_(x)+e)
140#define EFF_RAD_JAC_ELEM_(x,e) this%condensed_data_real(PHASE_REAL_LOC_(x)-1+e)
141#define NUM_CONC_JAC_ELEM_(x,e) this%condensed_data_real(PHASE_REAL_LOC_(x)-1+NUM_AERO_PHASE_JAC_ELEM_(x)+e)
142#define MASS_JAC_ELEM_(x,e) this%condensed_data_real(PHASE_REAL_LOC_(x)-1+2*NUM_AERO_PHASE_JAC_ELEM_(x)+e)
143#define MW_JAC_ELEM_(x,e) this%condensed_data_real(PHASE_REAL_LOC_(x)-1+3*NUM_AERO_PHASE_JAC_ELEM_(x)+e)
144
146
147 !> Generic test reaction data type
149 contains
150 !> Reaction initialization
151 procedure :: initialize
152 !> Finalize the reaction
155
156 !> Constructor for rxn_SIMPOL_phase_transfer_t
158 procedure :: constructor
159 end interface rxn_simpol_phase_transfer_t
160
161contains
162
163!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
164
165 !> Constructor for Phase transfer reaction
166 function constructor() result(new_obj)
167
168 !> A new reaction instance
169 type(rxn_simpol_phase_transfer_t), pointer :: new_obj
170
171 allocate(new_obj)
172 new_obj%rxn_phase = aero_rxn
173
174 end function constructor
175
176!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
177
178 !> Initialize the reaction data, validating component data and loading
179 !! any required information into the condensed data arrays for use during
180 !! solving
181 subroutine initialize(this, chem_spec_data, aero_rep, n_cells)
182
183 !> Reaction data
184 class(rxn_simpol_phase_transfer_t), intent(inout) :: this
185 !> Chemical species data
186 type(chem_spec_data_t), intent(in) :: chem_spec_data
187 !> Aerosol representations
188 type(aero_rep_data_ptr), pointer, intent(in) :: aero_rep(:)
189 !> Number of grid cells to solve simultaneously
190 integer(kind=i_kind), intent(in) :: n_cells
191
192 type(property_t), pointer :: spec_props, b_params
193 character(len=:), allocatable :: key_name, gas_spec_name, aero_spec_name
194 character(len=:), allocatable :: phase_name, act_name, error_msg
195 integer(kind=i_kind) :: i_spec, i_aero_rep, n_aero_ids, i_aero_id
196 integer(kind=i_kind) :: i_phase, n_aero_jac_elem, tmp_size
197 type(string_t), allocatable :: unique_spec_names(:), unique_act_names(:)
198 integer(kind=i_kind), allocatable :: phase_ids(:)
199 real(kind=dp) :: temp_real, n_star
200 logical :: has_act_coeff
201
202 ! Get the property set
203 if (.not. associated(this%property_set)) call die_msg(382913491, &
204 "Missing property set needed to initialize reaction")
205
206 ! Get the gas-phase species name
207 key_name = "gas-phase species"
208 call assert_msg(740333884, &
209 this%property_set%get_string(key_name, gas_spec_name), &
210 "Missing gas-phase species in SIMPOL.1 phase transfer reaction")
211
212 ! Get the aerosol phase name
213 key_name = "aerosol phase"
214 call assert_msg(325074932, &
215 this%property_set%get_string(key_name, phase_name), &
216 "Missing aerosol phase in SIMPOL.1 phase-transfer reaction")
217
218 ! Get the aerosol-phase species name
219 key_name = "aerosol-phase species"
220 call assert_msg(988456388, &
221 this%property_set%get_string(key_name, aero_spec_name), &
222 "Missing aerosol-phase species in SIMPOL.1 phase-transfer "// &
223 "reaction")
224
225 ! Set up a general error message
226 error_msg = " for SIMPOL.1 phase transfer of gas species '"// &
227 gas_spec_name//"' to aerosol-phase species '"// &
228 aero_spec_name//"' in phase '"//phase_name//"'"
229
230 ! Get the aerosol-phase activity coeffcient name
231 key_name = "aerosol-phase activity coefficient"
232 has_act_coeff = this%property_set%get_string(key_name, act_name)
233
234 ! Check for aerosol representations
235 call assert_msg(260518827, associated(aero_rep), &
236 "Missing aerosol representation"//error_msg)
237 call assert_msg(590304021, size(aero_rep).gt.0, &
238 "Missing aerosol representation"//error_msg)
239
240 ! Count the instances of this phase/species pair
241 n_aero_ids = 0
242 n_aero_jac_elem = 0
243 do i_aero_rep = 1, size(aero_rep)
244
245 ! Get the unique names in this aerosol representation for the
246 ! partitioning species
247 unique_spec_names = aero_rep(i_aero_rep)%val%unique_names( &
248 phase_name = phase_name, spec_name = aero_spec_name, &
249 phase_is_at_surface = .true.)
250
251 ! Skip aerosol representations that do not contain this phase
252 if (.not.allocated(unique_spec_names)) cycle
253
254 ! Add these instances to the list
255 n_aero_ids = n_aero_ids + size(unique_spec_names)
256
257 ! Get the number of Jacobian elements for calculations of mass, volume,
258 ! number, etc. for this partitioning into this phase
259 phase_ids = aero_rep(i_aero_rep)%val%phase_ids(phase_name, is_at_surface=.true.)
260 do i_phase = 1, size(phase_ids)
261 n_aero_jac_elem = n_aero_jac_elem + &
262 aero_rep(i_aero_rep)%val%num_jac_elem(phase_ids(i_phase))
263 end do
264
265 deallocate(unique_spec_names)
266
267 end do
268
269 call assert_msg(314000134, n_aero_ids.gt.0, &
270 "Aerosol species not found"//error_msg)
271
272 ! Allocate space in the condensed data arrays
273 allocate(this%condensed_data_int(num_int_prop_ + 2 + n_aero_ids * 13 + &
274 n_aero_jac_elem * 2))
275 allocate(this%condensed_data_real(num_real_prop_ + n_aero_jac_elem * 4))
276 this%condensed_data_int(:) = int(0, kind=i_kind)
277 this%condensed_data_real(:) = real(0.0, kind=dp)
278
279 ! Save space for the environment-dependent parameters
280 this%num_env_params = num_env_param_
281
282 ! Set the number of aerosol-species instances
283 num_aero_phase_ = n_aero_ids
284
285 ! Get the properties required of the aerosol species
286 call assert_msg(162662115, &
287 chem_spec_data%get_property_set(aero_spec_name, spec_props), &
288 "Missing properties"//error_msg)
289
290 ! Get the aerosol species molecular weight
291 key_name = "molecular weight [kg mol-1]"
292 call assert_msg(839930958, spec_props%get_real(key_name, mw_), &
293 "Missing property 'MW'"//error_msg)
294
295 ! Set the kg/m3 -> ppm conversion prefactor (multiply by T/P to get
296 ! conversion)
297 ! (ppm_x*Pa_air*m^3/K/kg_x) = Pa_air*m^3/mol_air/K * mol_x/kg_x *
298 ! 1.0e6ppm_x*mol_air/mol_x
299 conv_ = const%univ_gas_const / mw_ * 1.0e6
300
301 ! Set the ids of each aerosol-phase species instance
302 i_aero_id = 1
303 phase_int_loc_(i_aero_id) = num_int_prop_+12*num_aero_phase_+3
304 phase_real_loc_(i_aero_id) = num_real_prop_+1
305 do i_aero_rep = 1, size(aero_rep)
306
307 ! Get the unique names in this aerosol representation for the
308 ! partitioning species
309 unique_spec_names = aero_rep(i_aero_rep)%val%unique_names( &
310 phase_name = phase_name, spec_name = aero_spec_name, &
311 phase_is_at_surface = .true.)
312
313 ! Find the corresponding activity coefficients, if specified
314 if (has_act_coeff) then
315 unique_act_names = aero_rep(i_aero_rep)%val%unique_names( &
316 phase_name = phase_name, spec_name = act_name, &
317 phase_is_at_surface = .true.)
318 call assert_msg(236251734, size(unique_act_names).eq. &
319 size(unique_spec_names), &
320 "Mismatch of SIMPOL species and activity coeffs"// &
321 error_msg)
322 end if
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 ! Add the species concentration and activity coefficient ids to
327 ! the condensed data, and set the number of Jacobian elements for
328 ! the aerosol representations and the locations of the real data
329 do i_spec = 1, size(unique_spec_names)
330 num_aero_phase_jac_elem_(i_aero_id) = &
331 aero_rep(i_aero_rep)%val%num_jac_elem(phase_ids(i_spec))
332 aero_spec_(i_aero_id) = &
333 aero_rep(i_aero_rep)%val%spec_state_id( &
334 unique_spec_names(i_spec)%string)
335 if (has_act_coeff) then
336 aero_act_id_(i_aero_id) = &
337 aero_rep(i_aero_rep)%val%spec_state_id( &
338 unique_act_names(i_spec)%string)
339 else
340 aero_act_id_(i_aero_id) = -1
341 end if
342 aero_phase_id_(i_aero_id) = phase_ids(i_spec)
343 aero_rep_id_(i_aero_id) = i_aero_rep
344 i_aero_id = i_aero_id + 1
345 if (i_aero_id .le. num_aero_phase_) then
346 phase_int_loc_(i_aero_id) = phase_int_loc_(i_aero_id - 1) + 1 + &
347 2*num_aero_phase_jac_elem_(i_aero_id - 1)
348 phase_real_loc_(i_aero_id) = phase_real_loc_(i_aero_id - 1) + &
349 4*num_aero_phase_jac_elem_(i_aero_id - 1)
350 end if
351 end do
352
353 deallocate(unique_spec_names)
354
355 end do
356
357 ! Get the SIMPOL.1 parameters
358 key_name = "B"
359 call assert_msg(882881186, &
360 this%property_set%get_property_t(key_name, b_params), &
361 "Missing 'B' parameters"//error_msg)
362 call assert_msg(654885723, b_params%size().eq.4, &
363 "Incorrect number of 'B' parameters"//error_msg)
364 call b_params%iter_reset()
365 call assert_msg(694024883, b_params%get_real(val = b1_), &
366 "Got non-real 'B1' parameter"//error_msg)
367 call b_params%iter_next()
368 call assert_msg(231316411, b_params%get_real(val = b2_), &
369 "Got non-real 'B2' parameter"//error_msg)
370 call b_params%iter_next()
371 call assert_msg(126167907, b_params%get_real(val = b3_), &
372 "Got non-real 'B3' parameter"//error_msg)
373 call b_params%iter_next()
374 call assert_msg(573535753, b_params%get_real(val = b4_), &
375 "Got non-real 'B4' parameter"//error_msg)
376
377 ! Save the index of the gas-phase species in the state variable array
378 gas_spec_ = chem_spec_data%gas_state_id(gas_spec_name)
379
380 ! Make sure the species exists
381 call assert_msg(551477581, gas_spec_.gt.0, &
382 "Missing gas-phase species"//error_msg)
383
384 ! Get the required properties for the gas-phase species
385 call assert_msg(611221674, &
386 chem_spec_data%get_property_set(gas_spec_name, spec_props), &
387 "Missing properties for gas-phase species"//error_msg)
388
389 ! Get N* to calculate the mass accomodation coefficient. If it is not
390 ! present, set DELTA_H_ and DELTA_S_ to zero to indicate a mass
391 ! accomodation coefficient of 1.0
392 ! Mass accomodation equation is based on equations in:
393 ! Ervens, B., et al., 2003. "CAPRAM 2.4 (MODAC mechanism): An extended
394 ! and condensed tropospheric aqueous mechanism and its application."
395 ! J. Geophys. Res. 108, 4426. doi:10.1029/2002JD002202
396 key_name = "N star"
397 if (spec_props%get_real(key_name, n_star)) then
398 ! enthalpy change (kcal mol-1)
399 delta_h_ = real(- 10.0d0*(n_star-1.0d0) + &
400 7.53d0*(n_star**(2.0d0/3.0d0)-1.0d0) - 1.0d0, kind=dp)
401 ! entropy change (cal mol-1)
402 delta_s_ = real(- 32.0d0*(n_star-1.0d0) + &
403 9.21d0*(n_star**(2.0d0/3.0d0)-1.0d0) - 1.3d0, kind=dp)
404 ! Convert dH and dS to (J mol-1)
405 delta_h_ = real(delta_h_ * 4184.0d0, kind=dp)
406 delta_s_ = real(delta_s_ * 4.184d0, kind=dp)
407 else
408 delta_h_ = real(0.0, kind=dp)
409 delta_s_ = real(0.0, kind=dp)
410 end if
411
412 ! Get the diffusion coefficient (m^2/s)
413 key_name = "diffusion coeff [m2 s-1]"
414 call assert_msg(948176709, spec_props%get_real(key_name, diff_coeff_), &
415 "Missing diffusion coefficient"//error_msg)
416
417 ! Calculate the constant portion of c_rms [m /( K^(1/2) * s )]
418 key_name = "molecular weight [kg mol-1]"
419 call assert_msg(272813400, spec_props%get_real(key_name, temp_real), &
420 "Missing molecular weight"//error_msg)
421 pre_c_avg_ = sqrt(8.0*const%univ_gas_const/(const%pi*temp_real))
422
423 ! Check the sizes of the data arrays
424 tmp_size = phase_int_loc_(i_aero_id - 1) + 1 + &
425 2*num_aero_phase_jac_elem_(i_aero_id - 1) - 1
426 call assert_msg(625802519, size(this%condensed_data_int) .eq. tmp_size, &
427 "int array size mismatch"//error_msg)
428 tmp_size = phase_real_loc_(i_aero_id - 1) + &
429 4*num_aero_phase_jac_elem_(i_aero_id - 1) - 1
430 call assert_msg(391089510, size(this%condensed_data_real) .eq. tmp_size, &
431 "real array size mismatch"//error_msg)
432
433 end subroutine initialize
434
435!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
436
437 !> Finalize the reaction
438 subroutine finalize(this)
439
440 !> Reaction data
441 type(rxn_simpol_phase_transfer_t), intent(inout) :: this
442
443 if (associated(this%property_set)) &
444 deallocate(this%property_set)
445 if (allocated(this%condensed_data_real)) &
446 deallocate(this%condensed_data_real)
447 if (allocated(this%condensed_data_int)) &
448 deallocate(this%condensed_data_int)
449
450 end subroutine finalize
451
452!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
453
454 !> Finalize an array of reactions
455 subroutine finalize_array(this)
456
457 !> Array of reaction data
458 type(rxn_simpol_phase_transfer_t), intent(inout) :: this(:)
459
460 integer(kind=i_kind) :: i
461
462 do i = 1, size(this)
463 call finalize(this(i))
464 end do
465
466 end subroutine finalize_array
467
468!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
469
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:88
The rxn_SIMPOL_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 to aero_rep_data_t extending types.
Abstract reaction data type.
Definition rxn_data.F90:98
String type for building arrays of string of various size.
Definition util.F90:53