CAMP 1.0.0
Chemistry Across Multiple Phases
sub_model_PDFiTE.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_sub_model_PDFiTE module.
7!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
8!> \page camp_sub_model_PDFiTE CAMP: PD-FiTE Activity
9!!
10!! PD-FiTE activity calculates aerosol-phase species activities using
11!! Taylor series to describe partial derivatives of mean activity coefficients
12!! for ternary solutions, as described in Topping et al. (2009)
13!! \cite Topping2009 . The mean binary activity coefficients for ion pairs are
14!! calculated according to eq. 15 in \cite Topping2009 . The values are then
15!! made available to aqueous-phase reactions during solving.
16!!
17!! Input data for PDFiTE activity equations have the following format :
18!! \code{.json}
19!! { "camp-data" : [
20!! {
21!! "name" : "my PDFiTE activity",
22!! "type" : "SUB_MODEL_PDFITE",
23!! "gas-phase water" : "H2O",
24!! "aerosol-phase water" : "H2O_aq",
25!! "aerosol phase" : "my aero phase",
26!! "calculate for" : {
27!! "H-NO3" : {
28!! "interactions" : [
29!! {
30!! "ion pair" : "H-NO3",
31!! "min RH" : 0.0,
32!! "max RH" : 0.1,
33!! "B" : [ 0.925113 ]
34!! },
35!! {
36!! "ion pair" : "H-NO3",
37!! "min RH" : 0.1,
38!! "max RH" : 0.4,
39!! "B" : [ 0.12091, 13.497, -67.771, 144.01, -117.97 ]
40!! },
41!! {
42!! "ion pair" : "H-NO3",
43!! "min RH" : 0.4,
44!! "max RH" : 0.9,
45!! "B" : [ 1.3424, -0.8197, -0.52983, -0.37335 ]
46!! },
47!! {
48!! "ion pair" : "H-NO3",
49!! "min RH" : 0.9,
50!! "max RH" : 1.0,
51!! "B" : [ -0.3506505 ]
52!! },
53!! {
54!! "ion pair" : "NH4-NO3",
55!! "min RH" : 0.0,
56!! "max RH" : 0.1,
57!! "B" : [ -11.93308 ]
58!! },
59!! {
60!! "ion pair" : "NH4-NO3",
61!! "min RH" : 0.1,
62!! "max RH" : 0.99,
63!! "B" : [ -17.0372, 59.232, -86.312, 44.04 ]
64!! },
65!! {
66!! "ion pair" : "NH4-NO3",
67!! "min RH" : 0.99,
68!! "max RH" : 1.0,
69!! "B" : [ -0.2599432 ]
70!! }
71!! ]
72!! }
73!! ...
74!! }
75!! }
76!! ]}
77!! \endcode
78!! The key-value pair \b aerosol \b phase is required to specify the aerosol
79!! phase for which to calculate activity coefficients. The key-value pairs
80!! \b gas-phase \b water and \b aerosol-phase \b water must also be present
81!! and specify the names for the water species in each phase. The final
82!! required key-value pair is \b calculated \b for, which should contain a set
83!! of ion pairs that activity coefficients will be calculated for.
84!!
85!! The key names in this set must correspond to ion pairs that are present in
86!! the specified aerosol phase. The values must contain a key-value pair
87!! named \b interactions which includes an array of ion-pair interactions
88!! used to calculate equation 15 in \cite Topping2009 .
89!!
90!! Each element in the \b interactions array must include an \b ion pair
91!! that exists in the specified aerosol phase, a \b min \b RH and \b max \b RH
92!! that specify the bounds for which the fitted curve is valid, and an array
93!! of \b B values that specify the polynomial coefficients B0, B1, B2, ...
94!! as in equation 19 in \cite Topping2009 . At least one polynomial
95!! coefficient must be present.
96!!
97!! If at least one interaction with an ion pair is included, enough
98!! interactions with that ion pair must be included to cover the entire RH
99!! range (0.0--1.0). Interactions are assume to cover the range (minRH, maxRH],
100!! except for the lowest RH interaction, which covers th range [0.0, maxRH].
101!!
102!! When the interacting ion pair is the same as the ion-pair for which the
103!! mean binary activity coefficient is being calculated, the interaction
104!! parameters are used to calculate \f$ln(\gamma_A^0(\mathrm{RH}))\f$.
105!! Otherwise, the parameters are used to calculate
106!! \f$\frac{d\ln{\gamma_A(\mathrm{RH})}}{d(N_{B,M}N_{B,x})}\f$.
107!!
108!! For the above example, the following input data should be present:
109!! \code{.json}
110!! {
111!! "name" : "H2O",
112!! "type" : "CHEM_SPEC",
113!! "phase" : "GAS",
114!! },
115!! {
116!! "name" : "H2O_aq",
117!! "type" : "CHEM_SPEC",
118!! "phase" : "AEROSOL",
119!! },
120!! {
121!! "name" : "H_p",
122!! "type" : "CHEM_SPEC",
123!! "phase" : "AEROSOL",
124!! "charge" : 1,
125!! "molecular weight [kg mol-1]" : 1.008
126!! },
127!! {
128!! "name" : "NH4_p",
129!! "type" : "CHEM_SPEC",
130!! "phase" : "AEROSOL",
131!! "charge" : 1,
132!! "molecular weight" : 18.04
133!! },
134!! {
135!! "name" : "NO3_m",
136!! "type" : "CHEM_SPEC",
137!! "phase" : "AEROSOL",
138!! "charge" : -1
139!! "molecular weight [kg mol-1]" : 62.0049
140!! },
141!! {
142!! "name" : "NH4-NO3",
143!! "type" : "CHEM_SPEC",
144!! "tracer type" : "ION_PAIR",
145!! "ions" : {
146!! "NH4_p" : {},
147!! "NO3_m" : {}
148!! }
149!! },
150!! {
151!! "name" : "H-NO3",
152!! "type" : "CHEM_SPEC",
153!! "tracer type" : "ION_PAIR",
154!! "ions" : {
155!! "H_p" : {},
156!! "NO3_m" : {}
157!! }
158!! },
159!! {
160!! "name" : "my aero phase",
161!! "type" : "AERO_PHASE",
162!! "species" : ["H_p", "NO3_m", "NH4_p", "NH4-NO3", "H-NO3", "H2O_aq"]
163!! }
164!! \endcode
165!!
166!!
167
168!> The sub_model_PDFiTE_t type and associated functions.
170
171 use, intrinsic :: ieee_arithmetic, only: ieee_next_after
175 use camp_constants, only: const
177 use camp_property
179 use camp_util, only: i_kind, dp, to_string, &
181 die_msg, string_t, &
183
184 implicit none
185 private
186
187#define NUM_PHASE_ this%condensed_data_int(1)
188#define GAS_WATER_ID_ this%condensed_data_int(2)
189#define NUM_ION_PAIRS_ this%condensed_data_int(3)
190#define TOTAL_INT_PARAM_ this%condensed_data_int(4)
191#define TOTAL_FLOAT_PARAM_ this%condensed_data_int(5)
192#define NUM_INT_PROP_ 5
193#define NUM_REAL_PROP_ 0
194#define NUM_ENV_PARAM_ 1
195#define PHASE_ID_(x) this%condensed_data_int(NUM_INT_PROP_+x)
196#define PAIR_INT_PARAM_LOC_(x) this%condensed_data_int(NUM_INT_PROP_+NUM_PHASE_+x)
197#define PAIR_FLOAT_PARAM_LOC_(x) this%condensed_data_int(NUM_INT_PROP_+NUM_PHASE_+NUM_ION_PAIRS_+x)
198#define ION_PAIR_ACT_ID_(x) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x))
199#define NUM_CATION_(x) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+1)
200#define NUM_ANION_(x) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+2)
201#define CATION_ID_(x) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+3)
202#define ANION_ID_(x) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+4)
203#define NUM_INTER_(x) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+5)
204#define JAC_WATER_ID_(p,x) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+5+p)
205#define JAC_CATION_ID_(p,x,y) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+5+NUM_PHASE_+(p-1)*NUM_ION_PAIRS_+y)
206#define JAC_ANION_ID_(p,x,y) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+5+(1+NUM_ION_PAIRS_)*NUM_PHASE_+(p-1)*NUM_ION_PAIRS_+y)
207#define NUM_B_(x,y) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+5+(1+2*NUM_ION_PAIRS_)*NUM_PHASE_+y)
208#define INTER_SPEC_ID_(x,y) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+5+(1+2*NUM_ION_PAIRS_)*NUM_PHASE_+NUM_INTER_(x)+y)
209#define INTER_SPEC_LOC_(x,y) this%condensed_data_int(PAIR_INT_PARAM_LOC_(x)+5+(1+2*NUM_ION_PAIRS_)*NUM_PHASE_+2*(NUM_INTER_(x))+y)
210#define CATION_MW_(x) this%condensed_data_real(PAIR_FLOAT_PARAM_LOC_(x))
211#define ANION_MW_(x) this%condensed_data_real(PAIR_FLOAT_PARAM_LOC_(x)+1)
212#define CATION_N_(x) this%condensed_data_real(PAIR_FLOAT_PARAM_LOC_(x)+2)
213#define ANION_N_(x) this%condensed_data_real(PAIR_FLOAT_PARAM_LOC_(x)+3)
214#define MIN_RH_(x,y) this%condensed_data_real(INTER_SPEC_LOC_(x,y))
215#define MAX_RH_(x,y) this%condensed_data_real(INTER_SPEC_LOC_(x,y)+1)
216#define B_Z_(x,y,z) this%condensed_data_real(INTER_SPEC_LOC_(x,y)+1+z)
217
218 ! Update types (These must match values in sub_model_UNIFAC.c)
219 ! (none for now)
220
221 public :: sub_model_pdfite_t
222
223 !> Generic test reaction data type
224 type, extends(sub_model_data_t) :: sub_model_pdfite_t
225 contains
226 !> Reaction initialization
227 procedure :: initialize
228 !> Return a real number representing the priority of the sub-model
229 !! calculations. Low priority sub models may depend on the results
230 !! of higher priority sub models.
231 procedure :: priority
232 !> Finalize the reaction
234 end type sub_model_pdfite_t
235
236 !> Constructor for sub_model_PDFiTE_t
238 procedure :: constructor
239 end interface sub_model_pdfite_t
240
241contains
242
243!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
244
245 !> Constructor for activity reaction
246 function constructor() result(new_obj)
247
248 !> A new reaction instance
249 type(sub_model_pdfite_t), pointer :: new_obj
250
251 allocate(new_obj)
252
253 end function constructor
254
255!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
256
257 !> Initialize the reaction data, validating component data and loading
258 !! any required information into the condensed data arrays for use during
259 !! solving
260 subroutine initialize(this, aero_rep_set, aero_phase_set, chem_spec_data)
261
262 !> Reaction data
263 class(sub_model_pdfite_t), intent(inout) :: this
264 !> The set of aerosol representations
265 type(aero_rep_data_ptr), pointer, intent(in) :: aero_rep_set(:)
266 !> The set of aerosol phases
267 type(aero_phase_data_ptr), pointer, intent(in) :: aero_phase_set(:)
268 !> Chemical species data
269 type(chem_spec_data_t), intent(in) :: chem_spec_data
270
271 type(property_t), pointer :: spec_props, ion_pairs, ion_pair, &
272 sub_props, ions, interactions, interaction, poly_coeffs
273 character(len=:), allocatable :: key_name, spec_name, phase_name, &
274 string_val, inter_spec_name
275 integer(kind=i_kind) :: n_phase, n_ion_pair, n_int_param, n_float_param
276 integer(kind=i_kind) :: i_aero_rep, i_phase, i_ion_pair, i_ion, i_spec, &
277 i_poly_coeff, i_interaction, j_ion_pair, j_interaction
278 integer(kind=i_kind) :: qty, int_val, charge, total_charge, tracer_type
279 real(kind=dp) :: real_val, molecular_weight, min_rh, max_rh
280 type(string_t), allocatable :: unique_spec_names(:)
281 character(len=:), allocatable :: ion_pair_name, ion_name
282 type(string_t), allocatable :: ion_pair_names(:), temp_ion_pair_names(:)
283 integer(kind=i_kind), allocatable :: num_inter(:)
284 real(kind=dp), allocatable :: rh_range(:)
285
286 ! Get the reaction property set
287 if (.not. associated(this%property_set)) call die_msg(101529793, &
288 "Missing property set needed to initialize PDFiTE activity"// &
289 " reaction.")
290
291 ! Get the aerosol phase name
292 key_name = "aerosol phase"
293 call assert_msg(429861748, &
294 this%property_set%get_string(key_name, phase_name), &
295 "Missing aerosol phase in PDFiTE activity reaction.")
296
297 ! Count the instances of the aerosol phase
298 n_phase = 0
299 do i_aero_rep = 1, size(aero_rep_set)
300
301 ! Get the number of instances of the phase in this representation
302 n_phase = n_phase + &
303 aero_rep_set(i_aero_rep)%val%num_phase_instances(phase_name)
304
305 end do
306
307 call assert_msg(656403972, n_phase.gt.0, &
308 "Aerosol phase '"//phase_name//"' not present in any aerosol "// &
309 "representation for PDFiTE activity reaction.")
310
311 ! Get the ion pairs to calculate mean binary activity coefficients for
312 key_name = "calculate for"
313 call assert_msg(258251605, &
314 this%property_set%get_property_t(key_name, ion_pairs), &
315 "Missing ion pairs to calculate activity for in PDFiTE "// &
316 " activity reaction.")
317
318 ! Count the ion pairs to calculate activity coefficients for
319 n_ion_pair = ion_pairs%size()
320 call assert_msg(479529522, n_ion_pair .gt. 0, &
321 "Empty ion pair set in PDFiTE activity reaction.")
322
323 ! Allocate space for the ion pair names and array to check for their
324 ! interaction parameters
325 allocate(ion_pair_names(n_ion_pair))
326
327 ! Set the ion pair names
328 call ion_pairs%iter_reset()
329 do i_ion_pair = 1, n_ion_pair
330
331 ! Get the name of the ion_pair
332 call assert(497564051, ion_pairs%get_key(ion_pair_name))
333 ion_pair_names(i_ion_pair)%string = ion_pair_name
334
335 call ion_pairs%iter_next()
336 end do
337
338 ! Get the number of parameters required for each ion pair
339 ! Adding space for INT_PROPS and phase state ids
340 n_int_param = num_int_prop_ + n_phase
341 ! Adding space for REAL_PROPS
342 n_float_param = num_real_prop_
343 call ion_pairs%iter_reset()
344 do i_ion_pair = 1, n_ion_pair
345
346 ! Get the name of the ion_pair
347 call assert(680654801, ion_pairs%get_key(ion_pair_name))
348
349 ! Get the ion_pair properties
350 call assert_msg(287766741, ion_pairs%get_property_t(val=ion_pair), &
351 "Missing ion pair properties for '"//ion_pair_name// &
352 "' in PDFiTE activity reaction.")
353
354 ! Get the interactions
355 key_name = "interactions"
356 call assert_msg(883233319, &
357 ion_pair%get_property_t(key_name, interactions), &
358 "Missing interaction parameters for '"//ion_pair_name// &
359 "' in PDFiTE activity reaction.")
360
361 ! Loop through the interaction set and count the parameters for each
362 call interactions%iter_reset()
363 do i_interaction = 1, interactions%size()
364
365 ! Get the current interaction
366 call assert(347348920, interactions%get_property_t(val=interaction))
367
368 ! Get the name of the interacting ion pair for use in error messages
369 key_name = "ion pair"
370 call assert_msg(924945098, &
371 interaction%get_string(key_name, inter_spec_name), &
372 "Missing interaction species name for ion pair '"// &
373 ion_pair_name//"' in PD-FiTE activity reaction")
374
375 ! Get the set of polynomial coefficients
376 key_name = "B"
377 call assert_msg(157416082, &
378 interaction%get_property_t(key_name, poly_coeffs), &
379 "Missing 'B' array of polynomial coefficients for ion "// &
380 "pair '"//inter_spec_name//"' in activity coefficient "// &
381 "calculation for ion pair '"//ion_pair_name//"' in "// &
382 "PD-FiTE activity reaction.")
383
384 ! Check the number of polynomial coefficients
385 call assert_msg(320927773, poly_coeffs%size().gt.0, &
386 "Insufficient polynomial coefficients for PDFiTE activity "//&
387 "calculation for ion_pair '"//ion_pair_name//"'.")
388
389 ! Adding space for minRH, maxRH, B[]
390 n_float_param = n_float_param + 2 + poly_coeffs%size()
391
392 ! Adding space for size(B[]), interacting species id, location of
393 ! interaction parameters
394 n_int_param = n_int_param + 3
395
396 ! Make sure that this ion pair is in the list of ion pairs. If not,
397 ! add it to the list.
398 do j_ion_pair = 1, size(ion_pair_names)
399 if (inter_spec_name.eq.ion_pair_names(j_ion_pair)%string) exit
400 if (j_ion_pair.eq.size(ion_pair_names)) then
401 allocate(temp_ion_pair_names(j_ion_pair))
402 forall (i_spec=1:j_ion_pair)
403 temp_ion_pair_names(i_spec)%string = &
404 ion_pair_names(i_spec)%string
405 end forall
406 deallocate(ion_pair_names)
407 allocate(ion_pair_names(j_ion_pair+1))
408 forall (i_spec=1:j_ion_pair)
409 ion_pair_names(i_spec)%string = &
410 temp_ion_pair_names(i_spec)%string
411 end forall
412 deallocate(temp_ion_pair_names)
413 ion_pair_names(j_ion_pair+1)%string = inter_spec_name
414 end if
415 end do
416
417 call interactions%iter_next()
418 end do
419
420 call ion_pairs%iter_next()
421 end do
422
423 ! Update the number of ion pairs to include those from interactions that
424 ! activity coefficients are not calculated for
425 n_ion_pair = size(ion_pair_names)
426
427 ! Adding space for locations of int and float data for each ion pair,
428 ! ion pair activity coefficient state id, number of cations and anions,
429 ! state ids of cation and anion, the number of interactions,
430 ! Jacobian ids for each activity coeffcient per phase for contributions
431 ! from water and each other ion pair
432 n_int_param = n_int_param + (8+(1+2*n_ion_pair)*n_phase)*n_ion_pair
433 ! Adding space for MW of the cation and anion, and cation and anion
434 ! concentrations (for use during solving)
435 n_float_param = n_float_param + 4*n_ion_pair
436
437 ! Allocate space for an array used to make sure all ion pair interactions
438 ! are included and the RH range is covered for each
439 allocate(num_inter(n_ion_pair))
440 allocate(rh_range(n_ion_pair))
441
442 ! Allocate space in the condensed data arrays
443 allocate(this%condensed_data_int(n_int_param))
444 allocate(this%condensed_data_real(n_float_param))
445 this%condensed_data_int(:) = int(9999, kind=i_kind)
446 this%condensed_data_real(:) = real(9999.0, kind=dp)
447
448 ! Save space for the environment-dependent parameters
449 this%num_env_params = num_env_param_
450
451 ! Set some data dimensions
452 num_phase_ = n_phase
453 num_ion_pairs_ = n_ion_pair
454 total_int_param_ = n_int_param
455 total_float_param_ = n_float_param
456
457 ! Set the gas-phase water species
458 key_name = "gas-phase water"
459 call assert_msg(901506020, &
460 this%property_set%get_string(key_name, spec_name), &
461 "Missing gas-phase water species name in PDFiTE activity "// &
462 "reaction.")
463
464 gas_water_id_ = chem_spec_data%gas_state_id(spec_name)
465
466 call assert_msg(442608616, gas_water_id_ .gt. 0, &
467 "Cannot find gas-phase water species '"//spec_name//"' for "// &
468 "PDFiTE activity reaction.")
469
470 ! Set the aerosol-water species
471 key_name = "aerosol-phase water"
472 call assert_msg(214613153, &
473 this%property_set%get_string(key_name, spec_name), &
474 "Missing aerosol-phase water species name in PDFiTE activity "// &
475 "reaction.")
476
477 ! Make the PHASE_ID_(x) hold the state id of aerosol water in each
478 ! phase instance. Then the aerosol water id is 1, and the ion
479 ! ids will be relative to the water id in each phase.
480 i_phase = 1
481 do i_aero_rep = 1, size(aero_rep_set)
482 unique_spec_names = aero_rep_set(i_aero_rep)%val%unique_names( &
483 phase_name = phase_name, spec_name = spec_name)
484 if (.not.allocated(unique_spec_names)) cycle
485 do i_spec = 1, size(unique_spec_names)
486 phase_id_(i_phase) = aero_rep_set(i_aero_rep)%val%spec_state_id( &
487 unique_spec_names(i_spec)%string)
488 call assert(658622226, phase_id_(i_phase).gt.0)
489 i_phase = i_phase + 1
490 end do
491 deallocate(unique_spec_names)
492 end do
493 i_phase = i_phase - 1
494 call assert_msg(653357919, i_phase.eq.num_phase_, &
495 "Incorrect number of aerosol water instances in PDFiTE "// &
496 "activity reaction. Expected "//trim(to_string(num_phase_))// &
497 " but got "//trim(to_string(i_phase)))
498
499 ! Save the ion_pair parameters
500 n_int_param = num_int_prop_ + num_phase_ + 2*num_ion_pairs_
501 n_float_param = num_real_prop_
502 call ion_pairs%iter_reset()
503 do i_ion_pair = 1, n_ion_pair
504
505 ! Get the name of the ion_pair
506 ion_pair_name = ion_pair_names(i_ion_pair)%string
507
508 ! Set the location of this ion_pair's parameters in the condensed data
509 ! arrays.
510 pair_int_param_loc_(i_ion_pair) = n_int_param + 1
511 pair_float_param_loc_(i_ion_pair) = n_float_param + 1
512
513 ! Get the activity coefficient state ids for this ion_pair
514 i_phase = 1
515 do i_aero_rep = 1, size(aero_rep_set)
516 unique_spec_names = aero_rep_set(i_aero_rep)%val%unique_names( &
517 phase_name = phase_name, spec_name = ion_pair_name)
518 if (.not.allocated(unique_spec_names)) cycle
519 do i_spec = 1, size(unique_spec_names)
520 if (i_phase.eq.1) then
521 ion_pair_act_id_(i_ion_pair) = &
522 aero_rep_set(i_aero_rep)%val%spec_state_id( &
523 unique_spec_names(i_spec)%string) - &
524 phase_id_(i_phase)
525 else
526 call assert(142173386, ion_pair_act_id_(i_ion_pair).eq. &
527 aero_rep_set(i_aero_rep)%val%spec_state_id( &
528 unique_spec_names(i_spec)%string) - &
529 phase_id_(i_phase))
530 end if
531 i_phase = i_phase + 1
532 end do
533 deallocate(unique_spec_names)
534 end do
535 i_phase = i_phase - 1
536 call assert_msg(700406338, i_phase.eq.num_phase_, &
537 "Incorrect number of instances of ion pair '"// &
538 ion_pair_name//"' in PDFiTE activity reaction. Expected "// &
539 trim(to_string(num_phase_))//" but got "// &
540 trim(to_string(i_phase)))
541
542 ! Get the ion_pair species properties
543 call assert_msg(737691158, &
544 chem_spec_data%get_property_set(ion_pair_name, ion_pair), &
545 "Missing species properties for ion pair '"// &
546 ion_pair_name//"' in PDFiTE activity reaction.")
547
548 ! Make sure the specified ion pair is of the right tracer type
549 call assert(667033398, &
550 chem_spec_data%get_type(ion_pair_name, tracer_type))
551 call assert_msg(153203913, tracer_type.eq.chem_spec_activity_coeff, &
552 "Ion pair '"//ion_pair_name//"' must have a tracer type '"// &
553 "ION_PAIR' to participate in a PD-FiTE activity reaction.")
554
555 ! Get the number and id of the ions and their molecular weights
556 key_name = "ions"
557 call assert_msg(974596824, ion_pair%get_property_t(key_name, ions), &
558 "Mission ions for ion pair '"//ion_pair_name// &
559 "' in PDFiTE activity reaction.")
560
561 call assert_msg(852202160, ions%size().eq.2, &
562 "Invalid number of unique ions specified for ion pair '"// &
563 ion_pair_name//"' in for PDFiTE in activity reaction. "// &
564 "Expected 2 got "//trim(to_string(ions%size())))
565
566 ! TODO Consider moving some of the ion data validation to
567 ! camp_chem_spec_data
568 call ions%iter_reset()
569 total_charge = 0
570 do i_ion = 1, 2
571
572 ! Get the ion name
573 call assert(622753458, ions%get_key(ion_name))
574
575 ! Get the qty, if specified
576 qty = 1
577 if (ions%get_property_t(val=sub_props)) then
578 key_name = "qty"
579 if (sub_props%get_int(key_name, int_val)) qty = int_val
580 end if
581
582 ! Get the species properties
583 call assert_msg(619394685, &
584 chem_spec_data%get_property_set(ion_name, spec_props), &
585 "Missing species properties for ion '"//ion_name// &
586 "' in PDFiTE activity reaction.")
587
588 ! Get the molecular weight
589 key_name = "molecular weight [kg mol-1]"
590 call assert_msg(951085413, &
591 spec_props%get_real(key_name, molecular_weight), &
592 "Missing molecular weight for ion '"//ion_name// &
593 "' in PDFiTE activity reaction.")
594
595 ! Add the charge from this species
596 key_name = "charge"
597 call assert_msg(663798152, spec_props%get_int(key_name, charge), &
598 "Missing charge for ion '"//ion_name//"' in PDFiTE "// &
599 "activity reaction.")
600
601 if (charge.gt.0) then
602 num_cation_(i_ion_pair) = qty
603 cation_mw_(i_ion_pair) = molecular_weight
604 else if (charge.lt.0) then
605 num_anion_(i_ion_pair) = qty
606 anion_mw_(i_ion_pair) = molecular_weight
607 else
608 call die_msg(939555855, "Neutral species '"//ion_name// &
609 "' not allowed in PDFiTE activity reaction ion pair")
610 end if
611
612 ! Add contribution to total charge
613 total_charge = total_charge + qty * charge
614
615 ! Get the state ids for this species
616 i_phase = 1
617 do i_aero_rep = 1, size(aero_rep_set)
618 unique_spec_names = aero_rep_set(i_aero_rep)%val%unique_names( &
619 phase_name = phase_name, spec_name = ion_name)
620 if (.not.allocated(unique_spec_names)) cycle
621 do i_spec = 1, size(unique_spec_names)
622 if (charge.gt.0) then
623 if (i_phase.eq.1) then
624 cation_id_(i_ion_pair) = &
625 aero_rep_set(i_aero_rep)%val%spec_state_id( &
626 unique_spec_names(i_spec)%string) - &
627 phase_id_(i_phase)
628 else
629 call assert(425726370, cation_id_(i_ion_pair).eq. &
630 aero_rep_set(i_aero_rep)%val%spec_state_id( &
631 unique_spec_names(i_spec)%string) - &
632 phase_id_(i_phase))
633 end if
634 else
635 if (i_phase.eq.1) then
636 anion_id_(i_ion_pair) = &
637 aero_rep_set(i_aero_rep)%val%spec_state_id( &
638 unique_spec_names(i_spec)%string) - &
639 phase_id_(i_phase)
640 else
641 call assert(192466600, anion_id_(i_ion_pair).eq. &
642 aero_rep_set(i_aero_rep)%val%spec_state_id( &
643 unique_spec_names(i_spec)%string) - &
644 phase_id_(i_phase))
645 end if
646 end if
647 i_phase = i_phase + 1
648 end do
649 deallocate(unique_spec_names)
650 end do
651 i_phase = i_phase - 1
652 call assert_msg(759322632, i_phase.eq.num_phase_, &
653 "Incorrect number of instances of ion species '"// &
654 ion_name//"' in PDFiTE activity reaction. Expected "// &
655 trim(to_string(num_phase_))//" but got "// &
656 trim(to_string(i_phase)))
657
658 ! Get the next ion
659 call ions%iter_next()
660
661 end do
662
663 call assert_msg(415650051, total_charge.eq.0, &
664 "Charge imbalance for ion_pair '"//ion_pair_name// &
665 " in PDFiTE activity reaction. Total charge: "// &
666 trim(to_string(total_charge)))
667
668 n_float_param = n_float_param + 4
669 n_int_param = n_int_param + 6 + (1+2*num_ion_pairs_)*num_phase_
670
671 ! The first portion of ion_pair_names holds species for which
672 ! activity coefficients are calculated
673 if (i_ion_pair.le.ion_pairs%size()) then
674 call assert(362061689, ion_pairs%get_key(string_val))
675 call assert(444603372, string_val.eq.ion_pair_name)
676
677 ! Get the ion_pair properties
678 call assert(520886936, ion_pairs%get_property_t(val=ion_pair))
679
680 ! Get the interactions
681 key_name = "interactions"
682 call assert(216229321, ion_pair%get_property_t(key_name,interactions))
683
684 ! Set the number of interactions for this ion pair
685 num_inter_(i_ion_pair) = interactions%size()
686
687 ! Get the interaction parameters
688 num_inter(:) = 0
689 call interactions%iter_reset()
690 do i_interaction = 1, interactions%size()
691
692 ! Set the location of the interaction float parameters
693 inter_spec_loc_(i_ion_pair, i_interaction) = n_float_param + 1
694
695 ! Get the current interaction
696 call assert(105466070, interactions%get_property_t(val=interaction))
697
698 ! Get the name of the interacting species
699 key_name = "ion pair"
700 call assert_msg(220815620, &
701 interaction%get_string(key_name, inter_spec_name), &
702 "Missing interacting species name for ion pair '"// &
703 ion_pair_name//"' in PDFiTE activity reaction.")
704
705 ! Mark the species as having an interaction and get the id of the
706 ! ion_pair in the 'calculate for' list
707 do i_spec = 1, size(ion_pair_names)
708 if (ion_pair_names(i_spec)%string.eq.inter_spec_name) then
709 num_inter(i_spec) = num_inter(i_spec) + 1
710 inter_spec_id_(i_ion_pair,i_interaction) = i_spec
711 exit
712 end if
713 end do
714
715 ! Set the minimum RH value
716 key_name = "min RH"
717 call assert_msg(519674084, &
718 interaction%get_real(key_name, &
719 min_rh_(i_ion_pair, i_interaction)), &
720 "Missing minimum RH value for ion pair '"// &
721 ion_pair_name//"' interaction with '"//inter_spec_name// &
722 "' in PD-FiTE activity reaction.")
723 min_rh = min_rh_(i_ion_pair, i_interaction)
724 call assert_msg(294172408, &
725 min_rh.ge.real(0.0, kind=dp).and. &
726 min_rh.lt.real(1.0, kind=dp), &
727 "Invalid value for minimum RH for ion pair '"// &
728 ion_pair_name//"' interaction with '"//inter_spec_name// &
729 "' in PD-FiTE activity reaction: "//to_string(min_rh))
730
731 ! Set the maximum RH value
732 key_name = "max RH"
733 call assert_msg(649525712, &
734 interaction%get_real(key_name, &
735 max_rh_(i_ion_pair, i_interaction)), &
736 "Missing maximum RH value for ion pair '"// &
737 ion_pair_name//"' interaction with '"//inter_spec_name// &
738 "' in PD-FiTE activity reaction.")
739 max_rh = max_rh_(i_ion_pair, i_interaction)
740 call assert_msg(840423507, &
741 max_rh.gt.real(0.0, kind=dp).and. &
742 max_rh.le.real(1.0, kind=dp), &
743 "Invalid value for maximum RH for ion pair '"// &
744 ion_pair_name//"' interaction with '"//inter_spec_name// &
745 "' in PD-FiTE activity reaction: "//to_string(max_rh))
746
747 ! Get the number of B_z parameters
748 key_name = "B"
749 call assert(981216440, &
750 interaction%get_property_t(key_name, poly_coeffs))
751 num_b_(i_ion_pair, i_interaction) = poly_coeffs%size()
752
753 ! Get the B_z parameters
754 call poly_coeffs%iter_reset()
755 do i_poly_coeff = 1, poly_coeffs%size()
756 call assert_msg(988796931, poly_coeffs%get_real(val=real_val), &
757 "Invalid polynomial coefficient for ion_pair '"// &
758 ion_pair_name//"' interaction with '"//inter_spec_name// &
759 "' in PDFiTE activity reaction.")
760 b_z_(i_ion_pair, i_interaction, i_poly_coeff) = real_val
761 call poly_coeffs%iter_next()
762 end do
763
764 n_float_param = n_float_param + 2 + num_b_(i_ion_pair, i_interaction)
765 n_int_param = n_int_param + 3
766
767 call interactions%iter_next()
768 end do
769
770 ! Check that all the interactions were included at least once
771 do i_spec = 1, size(num_inter)
772 call warn_assert_msg(793223082, num_inter(i_spec).ge.1, &
773 "Missing interaction parameters between ion_pair '"// &
774 ion_pair_name//"' and '"//ion_pair_names(i_spec)%string// &
775 "' in PDFiTE activity interaction: "// &
776 trim(to_string(num_inter(i_spec)))//".")
777 end do
778
779 ! Make sure no interactions with the same ion pair overlap in
780 ! their RH ranges and that the entire RH range (0.0-1.0) is covered.
781 ! Treat ranges as R = (minRH,maxRH] to avoid overlap at boundaries
782 rh_range(:) = 0.0
783 do i_interaction = 1, num_inter_(i_ion_pair)
784
785 ! Check for RH range overlaps with other interactions with this
786 ! ion pair
787 do j_interaction = i_interaction+1, num_inter_(i_ion_pair)
788 if (inter_spec_id_(i_ion_pair, i_interaction).ne. &
789 inter_spec_id_(i_ion_pair, j_interaction)) cycle
790 call assert_msg(858420675, &
791 min_rh_(i_ion_pair, i_interaction).ge. &
792 max_rh_(i_ion_pair, j_interaction).or. &
793 max_rh_(i_ion_pair, i_interaction).le. &
794 min_rh_(i_ion_pair, j_interaction), &
795 "Overlapping RH range for interactions for ion pair '"// &
796 ion_pair_name//"' in PD-FiTE activity reaction.")
797 end do
798
799 ! Add the constribution from this interaction to RH coverage for
800 ! the interacting ion pair
801 rh_range(inter_spec_id_(i_ion_pair, i_interaction)) = &
802 rh_range(inter_spec_id_(i_ion_pair, i_interaction)) + &
803 (max_rh_(i_ion_pair, i_interaction) - &
804 min_rh_(i_ion_pair, i_interaction))
805 end do
806
807 ! Make sure the entire RH range is covered
808 do i_spec = 1, size(rh_range)
809 call assert_msg(370258071, &
810 (rh_range(i_spec).ge.ieee_next_after(1.0_dp, 0.0_dp) .and. &
811 rh_range(i_spec).le.ieee_next_after(1.0_dp, 2.0_dp)) .or. &
812 num_inter(i_spec).eq.0, &
813 "Incomplete RH coverage for interaction with ion pair '"// &
814 ion_pair_names(i_spec)%string//"' for '"//ion_pair_name// &
815 "' PD-FiTE activity coefficient calculation.")
816 end do
817
818 call ion_pairs%iter_next()
819
820 ! The last portion of ion_pair_names includes ion pairs that are
821 ! included in the interactions but for which acitivty coefficients
822 ! are not calculated
823 else
824
825 ! Set the number of interactions to zero to indicate not to calculate
826 ! activity coefficients for this ion pair
827 num_inter_(i_ion_pair) = 0
828
829 end if
830
831 end do
832
833 call assert(938415336, n_int_param.eq.total_int_param_)
834 call assert(433208931, n_float_param.eq.total_float_param_)
835
836 deallocate(ion_pair_names)
837 deallocate(num_inter)
838 deallocate(rh_range)
839
840 end subroutine initialize
841
842!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
843
844 !> Finalize the reaction
845 subroutine finalize(this)
846
847 !> Reaction data
848 type(sub_model_pdfite_t), intent(inout) :: this
849
850 if (associated(this%property_set)) &
851 deallocate(this%property_set)
852 if (allocated(this%condensed_data_real)) &
853 deallocate(this%condensed_data_real)
854 if (allocated(this%condensed_data_int)) &
855 deallocate(this%condensed_data_int)
856
857 end subroutine finalize
858
859!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
860
861 !> Finalize an array of reactions
862 subroutine finalize_array(this)
863
864 !> Array of reaction data
865 type(sub_model_pdfite_t), intent(inout) :: this(:)
866
867 integer(kind=i_kind) :: i
868
869 do i = 1, size(this)
870 call finalize(this(i))
871 end do
872
873 end subroutine finalize_array
874
875!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
876
877 !> Return a real number representing the priority of the sub model
878 !! calculations. Low priority sub models may use the results of higher
879 !! priority sub models. Lower numbers indicate higher priority.
880 !!
881 !! PD-FiTE calculations may depend on water concentrations, so must be
882 !! lower priority than the ZSR sub model.
883 function priority(this)
884
885 !> Sub model priority
886 real(kind=dp) :: priority
887 !> Sub model data
888 class(sub_model_pdfite_t), intent(in) :: this
889
890 priority = 2.0;
891
892 end function priority
893
894!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
895
896end module camp_sub_model_pdfite
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.
Return a real number representing the priority of the sub model calculations. Low priority sub models...
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.
integer(kind=i_kind), parameter, public chem_spec_activity_coeff
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 abstract sub_model_data_t structure and associated subroutines.
The sub_model_PDFiTE_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
subroutine warn_assert_msg(code, condition_ok, warning_msg)
Prints a warning message if condition_ok is false.
Definition util.F90:112
Pointer type for building arrays.
Pointer to aero_rep_data_t extending types.
Abstract sub-model data type.
String type for building arrays of string of various size.
Definition util.F90:53