CAMP 1.0.0
Chemistry Across Multiple Phases
rxn_wennberg_tunneling.F90
Go to the documentation of this file.
1! Copyright (C) 2021 Barcelona Supercomputing Center,
2! University of Illinois at Urbana-Champaign, and
3! National Center for Atmospheric Research
4! SPDX-License-Identifier: MIT
5
6!> \file
7!> The camp_rxn_wennberg_tunneling module.
8
9!> \page camp_rxn_wennberg_tunneling CAMP: Wennberg Tunneling Reaction
10!!
11!! Wennberg tunneling reaction rate constant equations are calculated as follows:
12!!
13!! \f[
14!! Ae^{(\frac{-B}{T})}e^{(\frac{C}{T^3})}
15!! \f]
16!!
17!! where \f$A\f$ is the pre-exponential factor
18!! (\f$(\mbox{#}\,\mbox{cm}^{-3})^{-(n-1)}\mbox{s}^{-1}\f$),
19!! and \f$B\f$ and \f$C\f$ are parameters that capture the temperature
20!! dependence as described in Wennberg et al. (2018) \cite Wennberg2018 .
21!!
22!! Input data for Wennberg tunneling equations has the following format:
23!! \code{.json}
24!! {
25!! "type" : "WENNBERG_TUNNELING",
26!! "A" : 123.45,
27!! "B" : 1200.0,
28!! "C" : 1.0e8,
29!! "time unit" : "MIN",
30!! "reactants" : {
31!! "spec1" : {},
32!! "spec2" : { "qty" : 2 },
33!! ...
34!! },
35!! "products" : {
36!! "spec3" : {},
37!! "spec4" : { "yield" : 0.65 },
38!! ...
39!! }
40!! }
41!! \endcode
42!! The key-value pairs \b reactants, and \b products are required. Reactants
43!! without a \b qty value are assumed to appear once in the reaction equation.
44!! Products without a specified \b yield are assumed to have a \b yield of
45!! 1.0.
46!!
47!! When \b A is not included, it is assumed to be 1.0, when \b B is not
48!! included, it is assumed to be 0.0 K, and when \b C is not included, it is
49!! assumed to be 0.0 \f$\mbox{K}^3\f$.
50!! The unit for time is assumed to be s, but inclusion of the optional
51!! key-value pair \b time \b unit = \b MIN can be used to indicate a rate
52!! with min as the time unit.
53
54!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
55
56!> The rxn_wennberg_tunneling_t type and associated functions.
58
61 use camp_constants, only: const
65 use camp_util, only: i_kind, dp, to_string, &
67
68 implicit none
69 private
70
71#define NUM_REACT_ this%condensed_data_int(1)
72#define NUM_PROD_ this%condensed_data_int(2)
73#define A_ this%condensed_data_real(1)
74#define B_ this%condensed_data_real(2)
75#define C_ this%condensed_data_real(3)
76#define CONV_ this%condensed_data_real(4)
77#define NUM_INT_PROP_ 2
78#define NUM_REAL_PROP_ 4
79#define NUM_ENV_PARAM_ 1
80#define REACT_(x) this%condensed_data_int(NUM_INT_PROP_ + x)
81#define PROD_(x) this%condensed_data_int(NUM_INT_PROP_ + NUM_REACT_ + x)
82#define DERIV_ID_(x) this%condensed_data_int(NUM_INT_PROP_ + NUM_REACT_ + NUM_PROD_ + x)
83#define JAC_ID_(x) this%condensed_data_int(NUM_INT_PROP_ + 2*(NUM_REACT_+NUM_PROD_) + x)
84#define YIELD_(x) this%condensed_data_real(NUM_REAL_PROP_ + x)
85
87
88 !> Generic test reaction data type
90 contains
91 !> Reaction initialization
92 procedure :: initialize
93 !> Finalize the reaction
94 final :: finalize
96
97 !> Constructor for rxn_wennberg_tunneling_t
99 procedure :: constructor
100 end interface rxn_wennberg_tunneling_t
101
102contains
103
104!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
105
106 !> Constructor for Wennberg tunneling reaction
107 function constructor() result(new_obj)
108
109 !> A new reaction instance
110 type(rxn_wennberg_tunneling_t), pointer :: new_obj
111
112 allocate(new_obj)
113 new_obj%rxn_phase = gas_rxn
114
115 end function constructor
116
117!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
118
119 !> Initialize the reaction data, validating component data and loading
120 !! any required information into the condensed data arrays for use during
121 !! solving
122 subroutine initialize(this, chem_spec_data, aero_rep, n_cells)
123
124 !> Reaction data
125 class(rxn_wennberg_tunneling_t), intent(inout) :: this
126 !> Chemical species data
127 type(chem_spec_data_t), intent(in) :: chem_spec_data
128 !> Aerosol representations
129 type(aero_rep_data_ptr), pointer, intent(in) :: aero_rep(:)
130 !> Number of grid cells being solved simultaneously
131 integer(kind=i_kind), intent(in) :: n_cells
132
133 type(property_t), pointer :: spec_props, reactants, products
134 character(len=:), allocatable :: key_name, spec_name, string_val
135 integer(kind=i_kind) :: i_spec, i_qty
136
137 integer(kind=i_kind) :: temp_int
138 real(kind=dp) :: temp_real
139
140 ! Get the species involved
141 if (.not. associated(this%property_set)) call die_msg(255324828, &
142 "Missing property set needed to initialize reaction")
143 key_name = "reactants"
144 call assert_msg(362119416, &
145 this%property_set%get_property_t(key_name, reactants), &
146 "Wennberg tunneling reaction is missing reactants")
147 key_name = "products"
148 call assert_msg(409429361, &
149 this%property_set%get_property_t(key_name, products), &
150 "Wennberg tunneling reaction is missing products")
151
152 ! Count the number of reactants (including those with a qty specified)
153 call reactants%iter_reset()
154 i_spec = 0
155 do while (reactants%get_key(spec_name))
156 ! Get properties included with this reactant in the reaction data
157 call assert(463909147, reactants%get_property_t(val=spec_props))
158 key_name = "qty"
159 if (spec_props%get_int(key_name, temp_int)) i_spec = i_spec+temp_int-1
160 call reactants%iter_next()
161 i_spec = i_spec + 1
162 end do
163
164 ! Allocate space in the condensed data arrays
165 allocate(this%condensed_data_int(num_int_prop_ + &
166 (i_spec + 2) * (i_spec + products%size())))
167 allocate(this%condensed_data_real(num_real_prop_ + products%size()))
168 this%condensed_data_int(:) = int(0, kind=i_kind)
169 this%condensed_data_real(:) = real(0.0, kind=dp)
170
171 ! Save space for the environment dependent parameters
172 this%num_env_params = num_env_param_
173
174 ! Save the size of the reactant and product arrays (for reactions where
175 ! these can vary)
176 num_react_ = i_spec
177 num_prod_ = products%size()
178
179 ! Set the #/cc -> ppm conversion prefactor
180 conv_ = const%avagadro / const%univ_gas_const * 10.0d0**(-12.0d0)
181
182 ! Get reaction parameters (it might be easiest to keep these at the
183 ! beginning of the condensed data array, so they can be accessed using
184 ! compliler flags)
185 key_name = "A"
186 if (.not. this%property_set%get_real(key_name, a_)) then
187 a_ = 1.0
188 end if
189 key_name = "time unit"
190 if (this%property_set%get_string(key_name, string_val)) then
191 if (trim(string_val).eq."MIN") then
192 a_ = a_ / 60.0
193 end if
194 endif
195 key_name = "B"
196 if (.not. this%property_set%get_real(key_name, b_)) then
197 b_ = 0.0
198 end if
199 key_name = "C"
200 if (.not. this%property_set%get_real(key_name, c_)) then
201 c_ = 0.0
202 end if
203
204 ! Get the indices and chemical properties for the reactants
205 call reactants%iter_reset()
206 i_spec = 1
207 do while (reactants%get_key(spec_name))
208
209 ! Save the index of this species in the state variable array
210 react_(i_spec) = chem_spec_data%gas_state_id(spec_name)
211
212 ! Make sure the species exists
213 call assert_msg(292299004, react_(i_spec).gt.0, &
214 "Missing Wennberg tunneling reactant: "//spec_name)
215
216 ! Get properties included with this reactant in the reaction data
217 call assert(687092598, reactants%get_property_t(val=spec_props))
218 key_name = "qty"
219 if (spec_props%get_int(key_name, temp_int)) then
220 do i_qty = 1, temp_int - 1
221 react_(i_spec + i_qty) = react_(i_spec)
222 end do
223 i_spec = i_spec + temp_int - 1
224 end if
225
226 call reactants%iter_next()
227 i_spec = i_spec + 1
228 end do
229
230 ! Get the indices and chemical properties for the products
231 call products%iter_reset()
232 i_spec = 1
233 do while (products%get_key(spec_name))
234
235 ! Save the index of this species in the state variable array
236 prod_(i_spec) = chem_spec_data%gas_state_id(spec_name)
237
238 ! Make sure the species exists
239 call assert_msg(681828291, prod_(i_spec).gt.0, &
240 "Missing Wennberg tunneling product: "//spec_name)
241
242 ! Get properties included with this product in the reaction data
243 call assert(794146636, products%get_property_t(val=spec_props))
244 key_name = "yield"
245 if (spec_props%get_real(key_name, temp_real)) then
246 yield_(i_spec) = temp_real
247 else
248 yield_(i_spec) = 1.0
249 end if
250
251 call products%iter_next()
252 i_spec = i_spec + 1
253 end do
254
255 end subroutine initialize
256
257!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
258
259 !> Finalize the reaction
260 elemental subroutine finalize(this)
261
262 !> Reaction data
263 type(rxn_wennberg_tunneling_t), intent(inout) :: this
264
265 if (associated(this%property_set)) &
266 deallocate(this%property_set)
267 if (allocated(this%condensed_data_real)) &
268 deallocate(this%condensed_data_real)
269 if (allocated(this%condensed_data_int)) &
270 deallocate(this%condensed_data_int)
271
272 end subroutine finalize
273
274!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
275
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_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
The rxn_wennberg_tunneling_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