CAMP 1.0.0
Chemistry Across Multiple Phases
rxn_CMAQ_H2O2.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_CMAQ_H2O2 module.
7
8!> \page camp_rxn_CMAQ_H2O2 CAMP: Special CMAQ Reaction for H2O2
9!!
10!! For the reactions:
11!!
12!! \f[\ce{
13!! HO2 + HO2 ->[k] H2O2
14!! }\f]
15!! \f[\ce{
16!! HO2 + HO2 + H2O ->[k] H2O2
17!! }\f]
18!!
19!! CMAQ rate constants are calculated as:
20!!
21!! \f[
22!! k=k_1+k_2[\mbox{M}]
23!! \f]
24!!
25!! where \f$k_1\f$ and \f$k_2\f$ are \ref camp_rxn_arrhenius "Arrhenius" rate
26!! constants with \f$D=300\f$ and \f$E=0\f$, and \f$[\mbox{M}]\f$ is the
27!! concentration of air \cite Gipson.
28!!
29!! Input data for CMAQ \f$\ce{H2O2}\f$ reactions have the following
30!! format:
31!! \code{.json}
32!! {
33!! "type" : "CMAQ_H2O2",
34!! "k1_A" : 5.6E-12,
35!! "k1_B" : -1.8,
36!! "k1_C" : 180.0,
37!! "k2_A" : 3.4E-12,
38!! "k2_B" : -1.6,
39!! "k2_C" : 104.1,
40!! "time unit" : "MIN",
41!! "reactants" : {
42!! "spec1" : {},
43!! "spec2" : { "qty" : 2 },
44!! ...
45!! },
46!! "products" : {
47!! "spec3" : {},
48!! "spec4" : { "yield" : 0.65 },
49!! ...
50!! }
51!! }
52!! \endcode
53!! The key-value pairs \b reactants, and \b products are required. Reactants
54!! without a \b qty value are assumed to appear once in the reaction equation.
55!! Products without a specified \b yield are assumed to have a \b yield of
56!! 1.0.
57!!
58!! The two sets of parameters beginning with \b k1_ and \b k2_ are the
59!! \ref camp_rxn_arrhenius "Arrhenius" parameters for the \f$k_1\f$ and
60!! \f$k_2\f$ rate constants, respectively. When not present, \b _A
61!! parameters are assumed to be 1.0, \b _B to be 0.0, and \b _C to be 0.0.
62!!
63!! The unit for time is assumed to be s, but inclusion of the optional
64!! key-value pair \b time \b unit = \b MIN can be used to indicate a rate
65!! with min as the time unit.
66
67!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
68
69!> The rxn_CMAQ_H2O2_t type and associated functions.
71
74 use camp_constants, only: const
78 use camp_util, only: i_kind, dp, to_string, &
80
81 implicit none
82 private
83
84#define NUM_REACT_ this%condensed_data_int(1)
85#define NUM_PROD_ this%condensed_data_int(2)
86#define k1_A_ this%condensed_data_real(1)
87#define k1_B_ this%condensed_data_real(2)
88#define k1_C_ this%condensed_data_real(3)
89#define k2_A_ this%condensed_data_real(4)
90#define k2_B_ this%condensed_data_real(5)
91#define k2_C_ this%condensed_data_real(6)
92#define CONV_ this%condensed_data_real(7)
93#define NUM_INT_PROP_ 2
94#define NUM_REAL_PROP_ 7
95#define NUM_ENV_PARAM_ 1
96#define REACT_(x) this%condensed_data_int(NUM_INT_PROP_ + x)
97#define PROD_(x) this%condensed_data_int(NUM_INT_PROP_ + NUM_REACT_ + x)
98#define DERIV_ID_(x) this%condensed_data_int(NUM_INT_PROP_ + NUM_REACT_ + NUM_PROD_ + x)
99#define JAC_ID_(x) this%condensed_data_int(NUM_INT_PROP_ + 2*(NUM_REACT_ + NUM_PROD_) + x)
100#define YIELD_(x) this%condensed_data_real(NUM_REAL_PROP_ + x)
101
102public :: rxn_cmaq_h2o2_t
103
104 !> Generic test reaction data type
105 type, extends(rxn_data_t) :: rxn_cmaq_h2o2_t
106 contains
107 !> Reaction initialization
108 procedure :: initialize
109 !> Finalize the reaction
110 final :: finalize
111 end type rxn_cmaq_h2o2_t
112
113 !> Constructor for rxn_CMAQ_H2O2_t
114 interface rxn_cmaq_h2o2_t
115 procedure :: constructor
116 end interface rxn_cmaq_h2o2_t
117
118contains
119
120!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
121
122 !> Constructor for CMAQ H2O2 reaction
123 function constructor() result(new_obj)
124
125 !> A new reaction instance
126 type(rxn_cmaq_h2o2_t), pointer :: new_obj
127
128 allocate(new_obj)
129 new_obj%rxn_phase = gas_rxn
130
131 end function constructor
132
133!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
134
135 !> Initialize the reaction data, validating component data and loading
136 !! any required information into the condensed data arrays for use during
137 !! solving
138 subroutine initialize(this, chem_spec_data, aero_rep, n_cells)
139
140 !> Reaction data
141 class(rxn_cmaq_h2o2_t), intent(inout) :: this
142 !> Chemical species data
143 type(chem_spec_data_t), intent(in) :: chem_spec_data
144 !> Aerosol representations
145 type(aero_rep_data_ptr), pointer, intent(in) :: aero_rep(:)
146 !> Number of grid cells to solve simultaneously
147 integer(kind=i_kind), intent(in) :: n_cells
148
149 type(property_t), pointer :: spec_props, reactants, products
150 character(len=:), allocatable :: key_name, spec_name, string_val
151 integer(kind=i_kind) :: i_spec, i_qty
152
153 integer(kind=i_kind) :: temp_int
154 real(kind=dp) :: temp_real
155
156 ! Get the species involved
157 if (.not. associated(this%property_set)) call die_msg(255324828, &
158 "Missing property set needed to initialize reaction")
159 key_name = "reactants"
160 call assert_msg(940945637, &
161 this%property_set%get_property_t(key_name, reactants), &
162 "CMAQ H2O2 reaction is missing reactants")
163 key_name = "products"
164 call assert_msg(435739232, &
165 this%property_set%get_property_t(key_name, products), &
166 "CMAQ H2O2 reaction is missing products")
167
168 ! Count the number of reactants (including those with a qty specified)
169 call reactants%iter_reset()
170 i_spec = 0
171 do while (reactants%get_key(spec_name))
172 ! Get properties included with this reactant in the reaction data
173 call assert(952475195, reactants%get_property_t(val=spec_props))
174 key_name = "qty"
175 if (spec_props%get_int(key_name, temp_int)) i_spec = i_spec+temp_int-1
176 call reactants%iter_next()
177 i_spec = i_spec + 1
178 end do
179
180 ! Allocate space in the condensed data arrays
181 ! Space in this example is allocated for two sets of inidices for the
182 ! reactants and products, one molecular property for each reactant,
183 ! yields for the products and three reaction parameters.
184 allocate(this%condensed_data_int(num_int_prop_ + &
185 (i_spec + 2) * (i_spec + products%size())))
186 allocate(this%condensed_data_real(num_real_prop_ + products%size()))
187 this%condensed_data_int(:) = int(0, kind=i_kind)
188 this%condensed_data_real(:) = real(0.0, kind=dp)
189
190 ! Save space for the environment dependent parameters
191 this%num_env_params = num_env_param_
192
193 ! Save the size of the reactant and product arrays (for reactions where
194 ! these can vary)
195 num_react_ = i_spec
196 num_prod_ = products%size()
197
198 ! Set the #/cc -> ppm conversion prefactor
199 conv_ = const%avagadro / const%univ_gas_const * 10.0d0**(-12.0d0)
200
201 ! Get reaction parameters (it might be easiest to keep these at the
202 ! beginning of the condensed data array, so they can be accessed using
203 ! compliler flags)
204 key_name = "k1_A"
205 if (.not. this%property_set%get_real(key_name, k1_a_)) then
206 k1_a_ = 1.0
207 end if
208 key_name = "k1_B"
209 if (.not. this%property_set%get_real(key_name, k1_b_)) then
210 k1_b_ = 0.0
211 end if
212 key_name = "k1_C"
213 if (.not. this%property_set%get_real(key_name, k1_c_)) then
214 k1_c_ = 0.0
215 end if
216 key_name = "k2_A"
217 if (.not. this%property_set%get_real(key_name, k2_a_)) then
218 k2_a_ = 1.0
219 end if
220 key_name = "k2_B"
221 if (.not. this%property_set%get_real(key_name, k2_b_)) then
222 k2_b_ = 0.0
223 end if
224 key_name = "k2_C"
225 if (.not. this%property_set%get_real(key_name, k2_c_)) then
226 k2_c_ = 0.0
227 end if
228 key_name = "time unit"
229 if (this%property_set%get_string(key_name, string_val)) then
230 if (trim(string_val).eq."MIN") then
231 k1_a_ = k1_a_ / 60.0
232 k2_a_ = k2_a_ / 60.0
233 end if
234 endif
235
236 ! Include the multiplication of [M]*k2 into k2_A_
237 k2_a_ = k2_a_ * real(1.0d6, kind=dp)
238
239 ! Get the indices and chemical properties for the reactants
240 call reactants%iter_reset()
241 i_spec = 1
242 do while (reactants%get_key(spec_name))
243
244 ! Save the index of this species in the state variable array
245 react_(i_spec) = chem_spec_data%gas_state_id(spec_name)
246
247 ! Make sure the species exists
248 call assert_msg(345360993, react_(i_spec).gt.0, &
249 "Missing CMAQ H2O2 reactant: "//spec_name)
250
251 ! Get properties included with this reactant in the reaction data
252 call assert(796763915, reactants%get_property_t(val=spec_props))
253 key_name = "qty"
254 if (spec_props%get_int(key_name, temp_int)) then
255 do i_qty = 1, temp_int - 1
256 react_(i_spec + i_qty) = react_(i_spec)
257 end do
258 i_spec = i_spec + temp_int - 1
259 end if
260
261 call reactants%iter_next()
262 i_spec = i_spec + 1
263 end do
264
265 ! Get the indices and chemical properties for the products
266 call products%iter_reset()
267 i_spec = 1
268 do while (products%get_key(spec_name))
269
270 ! Save the index of this species in the state variable array
271 prod_(i_spec) = chem_spec_data%gas_state_id(spec_name)
272
273 ! Make sure the species exists
274 call assert_msg(234948182, prod_(i_spec).gt.0, &
275 "Missing CMAQ H2O2 product: "//spec_name)
276
277 ! Get properties included with this product in the reaction data
278 call assert(267035567, products%get_property_t(val=spec_props))
279 key_name = "yield"
280 if (spec_props%get_real(key_name, temp_real)) then
281 yield_(i_spec) = temp_real
282 else
283 yield_(i_spec) = 1.0
284 end if
285
286 call products%iter_next()
287 i_spec = i_spec + 1
288 end do
289
290 end subroutine initialize
291
292!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
293
294 !> Finalize the reaction
295 elemental subroutine finalize(this)
296
297 !> Reaction data
298 type(rxn_cmaq_h2o2_t), intent(inout) :: this
299
300 if (associated(this%property_set)) &
301 deallocate(this%property_set)
302 if (allocated(this%condensed_data_real)) &
303 deallocate(this%condensed_data_real)
304 if (allocated(this%condensed_data_int)) &
305 deallocate(this%condensed_data_int)
306
307 end subroutine finalize
308
309!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
310
311end module camp_rxn_cmaq_h2o2
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_rep_data_t structure and associated subroutines.
The camp_state_t structure and associated subroutines.
Definition camp_state.F90:9
elemental subroutine finalize(this)
Finalize the state.
The chem_spec_data_t structure and associated subroutines.
type(chem_spec_data_t) function, pointer constructor(init_size)
Constructor for chem_spec_data_t.
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_CMAQ_H2O2_t type and associated functions.
The rxn_data_t structure and associated subroutines.
Definition rxn_data.F90:60
integer(kind=i_kind), parameter, public gas_rxn
Gas-phase reaction.
Definition rxn_data.F90:84
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.
Generic test reaction data type.
Abstract reaction data type.
Definition rxn_data.F90:98