CAMP 1.0.0
Chemistry Across Multiple Phases
rxn_wet_deposition.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_wet_deposition module.
7
8!> \page camp_rxn_wet_deposition CAMP: Wet Deposition
9!!
10!! Wet Deposition reactions take the form:
11!!
12!! \f[
13!! \mbox{X} \rightarrow
14!! \f]
15!!
16!! where \f$\ce{X}\f$ is a set of species in an aerosol phase
17!! undergoing wet deposition at a given rate.
18!!
19!! Wet deposition rate constants can be constant or set from an external
20!! module using the
21!! \c camp_rxn_wet_deposition::rxn_update_data_wet_deposition_t object.
22!! External modules should use the
23!! \c camp_rxn_wet_deposition::rxn_wet_deposition_t::get_property_set()
24!! function during initilialization to access any needed reaction parameters
25!! to identify certain wet deposition reactions.
26!! An \c camp_rxn_wet_deposition::update_data_wet_deposition_t object should be
27!! initialized for each wet deposition reaction. These objects can then be used
28!! during solving to update the wet deposition rate from an external module.
29!!
30!! Input data for wet deposition reactions have the following format :
31!! \code{.json}
32!! {
33!! "type" : "WET_DEPOSITION",
34!! "aerosol phase" : "my aero phase",
35!! "scaling factor" : 1.2,
36!! ...
37!! }
38!! \endcode
39!! The key-value pair \b aerosol \b phase is required and its value must be
40!! the name of the aerosol phase undergoing wet deposition. All species within
41!! the aerosol phase in all instances of the aerosol phase will be removed
42!! according the first-order loss rate constant. The \b scaling \b factor is
43!! optional, and can be used to set a constant scaling factor for the rate
44!! constant. When a \b scaling \b factor is not provided, it is assumed to be
45!! 1.0. All other data is optional and will be available to external modules
46!! during initialization. Rate constants are in units of \f$s^{-1}\f$, and
47!! must be set using a \c rxn_wet_deposition_update_data_t object.
48
49!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
50
51!> The rxn_wet_deposition_t type and associated functions.
53
56 use camp_constants, only: const
58 use camp_mpi
61 use camp_util, only: i_kind, dp, string_t, &
64
65 use iso_c_binding
66
67 implicit none
68 private
69
70#define RXN_ID_ this%condensed_data_int(1)
71#define NUM_SPEC_ this%condensed_data_int(2)
72#define SCALING_ this%condensed_data_real(1)
73#define NUM_INT_PROP_ 2
74#define NUM_REAL_PROP_ 1
75#define NUM_ENV_PARAM_ 2
76#define REACT_(s) this%condensed_data_int(NUM_INT_PROP_+s)
77#define DERIV_ID_(s) this%condensed_data_int(NUM_INT_PROP_+NUM_SPEC_+s)
78#define JAC_ID_(s) this%condensed_data_int(NUM_INT_PROP_+2*NUM_SPEC_+s))
79
81
82 !> Generic test reaction data type
83 type, extends(rxn_data_t) :: rxn_wet_deposition_t
84 contains
85 !> Reaction initialization
86 procedure :: initialize
87 !> Get the reaction property set
88 procedure :: get_property_set
89 !> Initialize update data
91 !> Finalize the reaction
94
95 !> Constructor for rxn_wet_deposition_t
97 procedure :: constructor
98 end interface rxn_wet_deposition_t
99
100 !> Wet Deposition rate update object
102 private
103 !> Flag indicating whether the update data as been allocated
104 logical :: is_malloced = .false.
105 !> Unique id for finding reactions during model initialization
106 integer(kind=i_kind) :: rxn_unique_id = 0
107 contains
108 !> Update the rate data
109 procedure :: set_rate => update_data_rate_set
110 !> Determine the pack size of the local update data
112 !> Pack the local update data to a binary
113 procedure :: internal_bin_pack
114 !> Unpack the local update data from a binary
116 !> Finalize the rate update data
119
120 !> Interface to c reaction functions
121 interface
122
123 !> Allocate space for a rate update
125 result(update_data) bind (c)
126 use iso_c_binding
127 !> Allocated update_data object
128 type(c_ptr) :: update_data
130
131 !> Set a new wet_deposition rate
133 rxn_unique_id, base_rate) bind (c)
134 use iso_c_binding
135 !> Update data
136 type(c_ptr), value :: update_data
137 !> Reaction id
138 integer(kind=c_int), value :: rxn_unique_id
139 !> New pre-scaling base wet_deposition rate
140 real(kind=c_double), value :: base_rate
142
143 !> Free an update rate data object
144 pure subroutine rxn_free_update_data(update_data) bind (c)
145 use iso_c_binding
146 !> Update data
147 type(c_ptr), value, intent(in) :: update_data
148 end subroutine rxn_free_update_data
149
150 end interface
151
152contains
153
154!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
155
156 !> Constructor for Wet Deposition reaction
157 function constructor() result(new_obj)
158
159 !> A new reaction instance
160 type(rxn_wet_deposition_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_rep, n_cells)
173
174 !> Reaction data
175 class(rxn_wet_deposition_t), intent(inout) :: this
176 !> Chemical species data
177 type(chem_spec_data_t), intent(in) :: chem_spec_data
178 !> Aerosol representations
179 type(aero_rep_data_ptr), pointer, intent(in) :: aero_rep(:)
180 !> Number of grid cells to solve simultaneously
181 integer(kind=i_kind), intent(in) :: n_cells
182
183 type(property_t), pointer :: spec_props
184 type(string_t), allocatable :: unique_names(:)
185 character(len=:), allocatable :: key_name, phase_name
186 integer(kind=i_kind) :: i_rep, i_spec, i_rep_spec, num_spec
187
188 integer(kind=i_kind) :: temp_int
189 real(kind=dp) :: temp_real
190
191 ! Get the reaction property set
192 call assert_msg(368664748, associated(this%property_set), &
193 "Missing property set needed to initialize reaction")
194
195 ! Get the aerosol phase name
196 key_name = "aerosol phase"
197 call assert_msg(133499444, &
198 this%property_set%get_string(key_name, phase_name), &
199 "Wet Deposition reaction is missing aerosol phase name")
200
201 ! Check for aerosol representations
202 call assert_msg(674938531, associated(aero_rep), &
203 "Missing aerosol representation for wet deposition reaction")
204 call assert_msg(731323851, size(aero_rep).gt.0, &
205 "Missing aerosol representation for wet deposition reaction")
206
207 ! Count the total number of species in the specified phase in each
208 ! aerosol representation
209 num_spec = 0
210 do i_rep = 1, size(aero_rep)
211 unique_names = aero_rep(i_rep)%val%unique_names( phase_name = &
212 phase_name )
213 num_spec = num_spec + size( unique_names )
214 end do
215 call assert_msg(332795980, num_spec.gt.0, &
216 "No species found for wet deposition aerosol phase "// &
217 phase_name)
218
219 ! Allocate space in the condensed data arrays
220 allocate(this%condensed_data_int(num_int_prop_+3*num_spec))
221 allocate(this%condensed_data_real(num_real_prop_))
222 this%condensed_data_int(:) = int(0, kind=i_kind)
223 this%condensed_data_real(:) = real(0.0, kind=dp)
224
225 ! Save space for the environment-dependent parameters
226 this%num_env_params = num_env_param_
227
228 ! Save the number of species
229 num_spec_ = num_spec
230
231 ! Get reaction parameters
232 key_name = "scaling factor"
233 if (.not. this%property_set%get_real(key_name, scaling_)) then
234 scaling_ = real(1.0, kind=dp)
235 end if
236
237 ! Save the indices of each species undergoing wet deposition
238 i_spec = 0
239 do i_rep = 1, size(aero_rep)
240 unique_names = aero_rep(i_rep)%val%unique_names( phase_name = &
241 phase_name )
242 do i_rep_spec = 1, size(unique_names)
243 i_spec = i_spec + 1
244 react_(i_spec) = aero_rep(i_rep)%val%spec_state_id( &
245 unique_names( i_rep_spec )%string )
246 call assert( 702159475, react_(i_spec) .gt. 0 )
247 end do
248 end do
249 call assert(312643342, i_spec .eq. num_spec_)
250
251 ! Initialize the reaction id
252 rxn_id_ = -1
253
254 end subroutine initialize
255
256!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
257
258 !> Get the reaction properties. (For use by external modules.)
259 function get_property_set(this) result(prop_set)
260
261 !> Reaction properties
262 type(property_t), pointer :: prop_set
263 !> Reaction data
264 class(rxn_wet_deposition_t), intent(in) :: this
265
266 prop_set => this%property_set
267
268 end function get_property_set
269
270!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
271
272 !> Finalize the reaction
273 subroutine finalize(this)
274
275 !> Reaction data
276 type(rxn_wet_deposition_t), intent(inout) :: this
277
278 if (associated(this%property_set)) &
279 deallocate(this%property_set)
280 if (allocated(this%condensed_data_real)) &
281 deallocate(this%condensed_data_real)
282 if (allocated(this%condensed_data_int)) &
283 deallocate(this%condensed_data_int)
284
285 end subroutine finalize
286
287!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
288
289 !> Finalize an array of reactions
290 subroutine finalize_array(this)
291
292 !> Array of reaction data
293 type(rxn_wet_deposition_t), intent(inout) :: this(:)
294
295 integer(kind=i_kind) :: i_rxn
296
297 do i_rxn = 1, size(this)
298 call finalize(this(i_rxn))
299 end do
300
301 end subroutine finalize_array
302
303!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
304
305 !> Set packed update data for wet_deposition rate constants
306 subroutine update_data_rate_set(this, base_rate)
307
308 !> Update data
309 class(rxn_update_data_wet_deposition_t), intent(inout) :: this
310 !> Updated pre-scaling wet_deposition rate
311 real(kind=dp), intent(in) :: base_rate
312
313 call rxn_wet_deposition_set_rate_update_data(this%get_data(), &
314 this%rxn_unique_id, base_rate)
315
316 end subroutine update_data_rate_set
317
318!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
319
320 !> Initialize update data
321 subroutine update_data_initialize(this, update_data, rxn_type)
322
323 use camp_rand, only : generate_int_id
324
325 !> The reaction to update
326 class(rxn_wet_deposition_t), intent(inout) :: this
327 !> Update data object
328 class(rxn_update_data_wet_deposition_t), intent(out) :: update_data
329 !> Reaction type id
330 integer(kind=i_kind), intent(in) :: rxn_type
331
332 ! If a reaction id has not yet been generated, do it now
333 if (rxn_id_.eq.-1) then
334 rxn_id_ = generate_int_id()
335 endif
336
337 update_data%rxn_unique_id = rxn_id_
338 update_data%rxn_type = int(rxn_type, kind=c_int)
339 update_data%update_data = rxn_wet_deposition_create_rate_update_data()
340 update_data%is_malloced = .true.
341
342 end subroutine update_data_initialize
343
344!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
345
346 !> Determine the size of a binary required to pack the reaction data
347 integer(kind=i_kind) function internal_pack_size(this, comm) &
348 result(pack_size)
349
350 !> Reaction update data
351 class(rxn_update_data_wet_deposition_t), intent(in) :: this
352 !> MPI communicator
353 integer, intent(in) :: comm
354
355 pack_size = &
356 camp_mpi_pack_size_logical(this%is_malloced, comm) + &
357 camp_mpi_pack_size_integer(this%rxn_unique_id, comm)
358
359 end function internal_pack_size
360
361!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
362
363 !> Pack the given value to the buffer, advancing position
364 subroutine internal_bin_pack(this, buffer, pos, comm)
365
366 !> Reaction update data
367 class(rxn_update_data_wet_deposition_t), intent(in) :: this
368 !> Memory buffer
369 character, intent(inout) :: buffer(:)
370 !> Current buffer position
371 integer, intent(inout) :: pos
372 !> MPI communicator
373 integer, intent(in) :: comm
374
375#ifdef CAMP_USE_MPI
376 integer :: prev_position
377
378 prev_position = pos
379 call camp_mpi_pack_logical(buffer, pos, this%is_malloced, comm)
380 call camp_mpi_pack_integer(buffer, pos, this%rxn_unique_id, comm)
381 call assert(865557010, &
382 pos - prev_position <= this%pack_size(comm))
383#endif
384
385 end subroutine internal_bin_pack
386
387!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
388
389 !> Unpack the given value from the buffer, advancing position
390 subroutine internal_bin_unpack(this, buffer, pos, comm)
391
392 !> Reaction update data
393 class(rxn_update_data_wet_deposition_t), intent(inout) :: this
394 !> Memory buffer
395 character, intent(inout) :: buffer(:)
396 !> Current buffer position
397 integer, intent(inout) :: pos
398 !> MPI communicator
399 integer, intent(in) :: comm
400
401#ifdef CAMP_USE_MPI
402 integer :: prev_position
403
404 prev_position = pos
405 call camp_mpi_unpack_logical(buffer, pos, this%is_malloced, comm)
406 call camp_mpi_unpack_integer(buffer, pos, this%rxn_unique_id, comm)
407 call assert(135713915, &
408 pos - prev_position <= this%pack_size(comm))
410#endif
411
412 end subroutine internal_bin_unpack
413
414!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
415
416 !> Finalize an update data object
417 subroutine update_data_finalize(this)
418
419 !> Update data object to free
420 type(rxn_update_data_wet_deposition_t), intent(inout) :: this
421
422 if (this%is_malloced) call rxn_free_update_data(this%update_data)
423
424 end subroutine update_data_finalize
425
426!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
427
428 !> Finalize an array of reactions
430
431 !> Array of update data objects
432 type(rxn_update_data_wet_deposition_t), intent(inout) :: this(:)
433
434 integer(kind=i_kind) :: i_rxn
435
436 do i_rxn = 1, size(this)
437 call update_data_finalize(this(i_rxn))
438 end do
439
440 end subroutine update_data_finalize_array
441
442!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
443
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 a list of unique names for each element on the camp_camp_state::camp_state_t::state_var array for...
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 aero_rxn
Aerosol-phase reaction.
Definition rxn_data.F90:88
The rxn_wet_deposition_t type and associated functions.
subroutine update_data_finalize(this)
Finalize an update data object.
subroutine update_data_initialize(this, update_data, rxn_type)
Initialize update data.
subroutine update_data_rate_set(this, base_rate)
Set packed update data for wet_deposition rate constants.
subroutine update_data_finalize_array(this)
Finalize an array of reactions.
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