CAMP 1.0.0
Chemistry Across Multiple Phases
rxn_first_order_loss.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_first_order_loss module.
7!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
8!> \page camp_rxn_first_order_loss CAMP: First-Order Loss
9!!
10!! First-Order Loss reactions take the form:
11!!
12!! \f[
13!! \mbox{X} \rightarrow
14!! \f]
15!!
16!! where \f$\ce{X}\f$ is the species being lost.
17!!
18!! First-Order loss rate constants can be constant or set from an external
19!! module using the
20!! \c camp_rxn_first_order_loss::rxn_update_data_first_order_loss_t object.
21!! External modules can use the
22!! \c camp_rxn_first_order_loss::rxn_first_order_loss_t::get_property_set()
23!! function during initilialization to access any needed reaction parameters
24!! to identify certain first-order loss reactions.
25!! An \c camp_rxn_first_order_loss::update_data_first_order_loss_t object
26!! should be initialized for each reaction. These objects can then
27!! be used during solving to update the first order loss rate from an
28!! external module.
29!!
30!! Input data for first-order loss reactions have the following format :
31!! \code{.json}
32!! {
33!! "type" : "FIRST_ORDER_LOSS",
34!! "species" : "species_name",
35!! "scaling factor" : 1.2,
36!! ...
37!! }
38!! \endcode
39!! The key-value pairs \b species is required and its value must be the name
40!! of the species being removed by the reaction. The \b scaling \b factor is
41!! optional, and
42!! can be used to set a constant scaling factor for the rate constant. When a
43!! \b scaling \b factor is not provided, it is assumed to be 1.0. All other
44!! data is optional and will be available to external modules during
45!! initialization. Rate constants are in units of \f$s^{-1}\f$, and must be
46!! set using a \c rxn_first_order_loss_update_data_t object.
47!!
48!!
49
50!> The rxn_first_order_loss_t type and associated functions.
52
55 use camp_constants, only: const
57 use camp_mpi
60 use camp_util, only: i_kind, dp, to_string, &
62
63 use iso_c_binding
64
65 implicit none
66 private
67
68#define RXN_ID_ this%condensed_data_int(1)
69#define REACT_ this%condensed_data_int(2)
70#define DERIV_ID_ this%condensed_data_int(3)
71#define JAC_ID_ this%condensed_data_int(4)
72#define SCALING_ this%condensed_data_real(1)
73#define NUM_INT_PROP_ 4
74#define NUM_REAL_PROP_ 1
75#define NUM_ENV_PARAM_ 2
76
78
79 !> Generic test reaction data type
80 type, extends(rxn_data_t) :: rxn_first_order_loss_t
81 contains
82 !> Reaction initialization
83 procedure :: initialize
84 !> Get the reaction property set
85 procedure :: get_property_set
86 !> Initialize update data
88 !> Finalize the reaction
91
92 !> Constructor for rxn_first_order_loss_t
94 procedure :: constructor
95 end interface rxn_first_order_loss_t
96
97 !> First-Order Loss rate update object
99 private
100 !> Flag indicating whether the update data as been allocated
101 logical :: is_malloced = .false.
102 !> Unique id for finding reactions during model initialization
103 integer(kind=i_kind) :: rxn_unique_id = 0
104 contains
105 !> Update the rate data
106 procedure :: set_rate => update_data_rate_set
107 !> Determine the pack size of the local update data
109 !> Pack the local update data to a binary
110 procedure :: internal_bin_pack
111 !> Unpack the local update data from a binary
113 !> Finalize the rate update data
116
117 !> Interface to c reaction functions
118 interface
119
120 !> Allocate space for a rate update
122 result(update_data) bind (c)
123 use iso_c_binding
124 !> Allocated update_data object
125 type(c_ptr) :: update_data
127
128 !> Set a new first_order_loss rate
130 rxn_unique_id, base_rate) bind (c)
131 use iso_c_binding
132 !> Update data
133 type(c_ptr), value :: update_data
134 !> Reaction id
135 integer(kind=c_int), value :: rxn_unique_id
136 !> New pre-scaling base first_order_loss rate
137 real(kind=c_double), value :: base_rate
139
140 !> Free an update rate data object
141 pure subroutine rxn_free_update_data(update_data) bind (c)
142 use iso_c_binding
143 !> Update data
144 type(c_ptr), value, intent(in) :: update_data
145 end subroutine rxn_free_update_data
146
147 end interface
148
149contains
150
151!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
152
153 !> Constructor for First-Order Loss reaction
154 function constructor() result(new_obj)
155
156 !> A new reaction instance
157 type(rxn_first_order_loss_t), pointer :: new_obj
158
159 allocate(new_obj)
160 new_obj%rxn_phase = gas_rxn
161
162 end function constructor
163
164!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
165
166 !> Initialize the reaction data, validating component data and loading
167 !! any required information into the condensed data arrays for use during
168 !! solving
169 subroutine initialize(this, chem_spec_data, aero_rep, n_cells)
170
171 !> Reaction data
172 class(rxn_first_order_loss_t), intent(inout) :: this
173 !> Chemical species data
174 type(chem_spec_data_t), intent(in) :: chem_spec_data
175 !> Aerosol representations
176 type(aero_rep_data_ptr), pointer, intent(in) :: aero_rep(:)
177 !> Number of grid cells to solve simultaneously
178 integer(kind=i_kind), intent(in) :: n_cells
179
180 type(property_t), pointer :: spec_props
181 character(len=:), allocatable :: key_name, spec_name
182 integer(kind=i_kind) :: i_spec, i_qty
183
184 integer(kind=i_kind) :: temp_int
185 real(kind=dp) :: temp_real
186
187 ! Get the species involved
188 call assert_msg(128411383, associated(this%property_set), &
189 "Missing property set needed to initialize reaction")
190 key_name = "species"
191 call assert_msg(164644065, &
192 this%property_set%get_string(key_name, spec_name), &
193 "First-Order Loss reaction is missing species name")
194
195 ! Allocate space in the condensed data arrays
196 allocate(this%condensed_data_int(num_int_prop_))
197 allocate(this%condensed_data_real(num_real_prop_))
198 this%condensed_data_int(:) = int(0, kind=i_kind)
199 this%condensed_data_real(:) = real(0.0, kind=dp)
200
201 ! Save space for the environment-dependent parameters
202 this%num_env_params = num_env_param_
203
204 ! Get reaction parameters
205 key_name = "scaling factor"
206 if (.not. this%property_set%get_real(key_name, scaling_)) then
207 scaling_ = real(1.0, kind=dp)
208 end if
209
210 ! Save the index of this species in the state variable array
211 react_ = chem_spec_data%gas_state_id(spec_name)
212
213 ! Make sure the species exists
214 call assert_msg(331442196, react_.gt.0, &
215 "Missing first-order loss species: "//spec_name)
216
217 ! Initialize the reaction id
218 rxn_id_ = -1
219
220 end subroutine initialize
221
222!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
223
224 !> Get the reaction properties. (For use by external modules.)
225 function get_property_set(this) result(prop_set)
226
227 !> Reaction properties
228 type(property_t), pointer :: prop_set
229 !> Reaction data
230 class(rxn_first_order_loss_t), intent(in) :: this
231
232 prop_set => this%property_set
233
234 end function get_property_set
235
236!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
237
238 !> Finalize the reaction
239 subroutine finalize(this)
240
241 !> Reaction data
242 type(rxn_first_order_loss_t), intent(inout) :: this
243
244 if (associated(this%property_set)) &
245 deallocate(this%property_set)
246 if (allocated(this%condensed_data_real)) &
247 deallocate(this%condensed_data_real)
248 if (allocated(this%condensed_data_int)) &
249 deallocate(this%condensed_data_int)
250
251 end subroutine finalize
252
253!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
254
255 !> Finalize an array of reactions
256 subroutine finalize_array(this)
257
258 !> Array of reaction data
259 type(rxn_first_order_loss_t), intent(inout) :: this(:)
260
261 integer(kind=i_kind) :: i
262
263 do i = 1, size(this)
264 call finalize(this(i))
265 end do
266
267 end subroutine finalize_array
268
269!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
270
271 !> Set packed update data for first_order_loss rate constants
272 subroutine update_data_rate_set(this, base_rate)
273
274 !> Update data
275 class(rxn_update_data_first_order_loss_t), intent(inout) :: this
276 !> Updated pre-scaling first_order_loss rate
277 real(kind=dp), intent(in) :: base_rate
278
279 call rxn_first_order_loss_set_rate_update_data(this%get_data(), &
280 this%rxn_unique_id, base_rate)
281
282 end subroutine update_data_rate_set
283
284!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
285
286 !> Initialize update data
287 subroutine update_data_initialize(this, update_data, rxn_type)
288
289 use camp_rand, only : generate_int_id
290
291 !> The reaction to update
292 class(rxn_first_order_loss_t), intent(inout) :: this
293 !> Update data object
294 class(rxn_update_data_first_order_loss_t), intent(out) :: update_data
295 !> Reaction type id
296 integer(kind=i_kind), intent(in) :: rxn_type
297
298 ! If a reaction id has not yet been generated, do it now
299 if (rxn_id_.eq.-1) then
300 rxn_id_ = generate_int_id()
301 endif
302
303 update_data%rxn_unique_id = rxn_id_
304 update_data%rxn_type = int(rxn_type, kind=c_int)
305 update_data%update_data = rxn_first_order_loss_create_rate_update_data()
306 update_data%is_malloced = .true.
307
308 end subroutine update_data_initialize
309
310!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
311
312 !> Determine the size of a binary required to pack the reaction data
313 integer(kind=i_kind) function internal_pack_size(this, comm) &
314 result(pack_size)
315
316 !> Reaction update data
317 class(rxn_update_data_first_order_loss_t), intent(in) :: this
318 !> MPI communicator
319 integer, intent(in) :: comm
320
321 pack_size = &
322 camp_mpi_pack_size_logical(this%is_malloced, comm) + &
323 camp_mpi_pack_size_integer(this%rxn_unique_id, comm)
324
325 end function internal_pack_size
326
327!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
328
329 !> Pack the given value to the buffer, advancing position
330 subroutine internal_bin_pack(this, buffer, pos, comm)
331
332 !> Reaction update data
333 class(rxn_update_data_first_order_loss_t), intent(in) :: this
334 !> Memory buffer
335 character, intent(inout) :: buffer(:)
336 !> Current buffer position
337 integer, intent(inout) :: pos
338 !> MPI communicator
339 integer, intent(in) :: comm
340
341#ifdef CAMP_USE_MPI
342 integer :: prev_position
343
344 prev_position = pos
345 call camp_mpi_pack_logical(buffer, pos, this%is_malloced, comm)
346 call camp_mpi_pack_integer(buffer, pos, this%rxn_unique_id, comm)
347 call assert(373785697, &
348 pos - prev_position <= this%pack_size(comm))
349#endif
350
351 end subroutine internal_bin_pack
352
353!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
354
355 !> Unpack the given value from the buffer, advancing position
356 subroutine internal_bin_unpack(this, buffer, pos, comm)
357
358 !> Reaction update data
359 class(rxn_update_data_first_order_loss_t), intent(inout) :: this
360 !> Memory buffer
361 character, intent(inout) :: buffer(:)
362 !> Current buffer position
363 integer, intent(inout) :: pos
364 !> MPI communicator
365 integer, intent(in) :: comm
366
367#ifdef CAMP_USE_MPI
368 integer :: prev_position
369
370 prev_position = pos
371 call camp_mpi_unpack_logical(buffer, pos, this%is_malloced, comm)
372 call camp_mpi_unpack_integer(buffer, pos, this%rxn_unique_id, comm)
373 call assert(368521390, &
374 pos - prev_position <= this%pack_size(comm))
376#endif
377
378 end subroutine internal_bin_unpack
379
380!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
381
382 !> Finalize an update data object
383 subroutine update_data_finalize(this)
384
385 !> Update data object to free
386 type(rxn_update_data_first_order_loss_t), intent(inout) :: this
387
388 if (this%is_malloced) call rxn_free_update_data(this%update_data)
389
390 end subroutine update_data_finalize
391
392!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
393
394 !> Finalize an array of update data objects
396
397 !> Array of update data objects to free
398 type(rxn_update_data_first_order_loss_t), intent(inout) :: this(:)
399
400 integer(kind=i_kind) :: i
401
402 do i = 1, size(this)
403 call update_data_finalize(this(i))
404 end do
405
406 end subroutine update_data_finalize_array
407
408!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
409
Initialize the aerosol representation data, validating component data and loading any required inform...
Extending-type binary pack function (Internal use only)
Extending-type binary unpack function (Internal use only)
Extending-type binary pack size (internal use only)
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.
integer(kind=i_kind) function pack_size(this, comm)
Determine the size of a binary required to pack the aerosol representation data.
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.
subroutine finalize_array(this_array)
Finalize the chemical species data.
type(chem_spec_data_t) function, pointer constructor(init_size)
Constructor for chem_spec_data_t.
logical function get_property_set(this, spec_name, property_set)
Get a species property set. Returns true if the species is found, or false otherwise.
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
Wrapper functions for MPI.
Definition mpi.F90:13
subroutine camp_mpi_pack_logical(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition mpi.F90:792
subroutine camp_mpi_unpack_integer(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition mpi.F90:1023
subroutine camp_mpi_unpack_logical(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition mpi.F90:1131
integer function camp_mpi_pack_size_logical(val, comm)
Determines the number of bytes required to pack the given value.
Definition mpi.F90:484
subroutine camp_mpi_pack_integer(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition mpi.F90:691
integer function camp_mpi_pack_size_integer(val, comm)
Determines the number of bytes required to pack the given value.
Definition mpi.F90:398
The property_t structure and associated subroutines.
Definition property.F90:9
Random number generators.
Definition rand.F90:9
integer(kind=i_kind) function generate_int_id()
Generate an integer id Ids will be sequential, and can only be generated by the primary process.
Definition rand.F90:435
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_first_order_loss_t type and associated functions.
subroutine update_data_finalize_array(this)
Finalize an array of update data objects.
subroutine update_data_initialize(this, update_data, rxn_type)
Initialize update data.
subroutine update_data_finalize(this)
Finalize an update data object.
subroutine update_data_rate_set(this, base_rate)
Set packed update data for first_order_loss rate constants.
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