CAMP 1.0.0
Chemistry Across Multiple Phases
chem_spec_data.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_chem_spec_data module.
7
8!> \page camp_species CAMP: Chemical Species
9!!
10!! Chemical species in the \ref index "camp-chem" module are gas- or
11!! aerosol-phase species that can participate in chemical or phase-transfer
12!! reactions in the \ref camp_mechanism "mechanism". Each species must
13!! have a unique name (the same name cannot be used for a gas-phase
14!! species and an aerosol-phase species). Chemical species may be present
15!! either in the gas-phase or any number of \ref camp_aero_phase
16!! "aerosol phases". For example, an aerosol-phase chemical species may be
17!! present in an "organic" and an "aqueous" \ref camp_aero_phase
18!! "aerosol phase".
19!!
20!! Chemical species data include physical constants and species-specific
21!! model parameters that are used during initialization to assemble reaction
22!! and sub-model data for use during solving. Note that chemical species data
23!! are **only** available during initialization, and when using MPI are not
24!! passed to child nodes. The primary node will, however, have access to data
25!! in the \c camp_camp_core::camp_core_t::chem_spec_data object for outputing
26!! model data (e.g., species names).
27!!
28!! The input format for chemical species can be found \ref
29!! input_format_species "here".
30
31!> The chem_spec_data_t structure and associated subroutines.
33
34#ifdef CAMP_USE_JSON
35 use json_module
36#endif
37 use camp_constants, only : dp, i_kind
40
41 use iso_c_binding
42
43 implicit none
44 private
45
46 public :: chem_spec_data_t
47
48 !> State variable types (Must match values in camp_solver.c)
49 integer(kind=i_kind), parameter, public :: chem_spec_unknown_type = 0
50 integer(kind=i_kind), parameter, public :: chem_spec_variable = 1
51 integer(kind=i_kind), parameter, public :: chem_spec_constant = 2
52 integer(kind=i_kind), parameter, public :: chem_spec_pssa = 3
53 integer(kind=i_kind), parameter, public :: chem_spec_activity_coeff = 4
54
55 !> Species phase
56 integer(kind=i_kind), parameter, public :: chem_spec_unknown_phase = 0
57 integer(kind=i_kind), parameter, public :: chem_spec_gas_phase = 1
58 integer(kind=i_kind), parameter, public :: chem_spec_aero_phase = 2
59
60 !> Reallocation increment
61 integer(kind=i_kind), parameter :: realloc_inc = 50
62 !> Default absolute integration tolerance
63 real(kind=dp), parameter :: default_abs_tol = 1.0e-14
64
65 !> Chemical species data
66 !!
67 !! Time-invariant data related to a \ref camp_species "chemical species"
69 private
70 !> Number of species
71 integer(kind=i_kind) :: num_spec = 0
72 !> Species name
73 type(string_t), pointer :: spec_name(:) => null()
74 !> Species type
75 integer(kind=i_kind), pointer :: spec_type(:) => null()
76 !> Species phase
77 integer(kind=i_kind), pointer :: spec_phase(:) => null()
78 !> Species property set
79 type(property_t), pointer :: property_set(:) => null()
80 contains
81 !> Load species from an input file
82 procedure :: load
83 !> Initialize the species set
84 procedure :: initialize
85 !> Get the number of species with specified conditions
86 procedure :: size => get_size
87 !> Check if a species name is in the set of chemical species
88 procedure :: exists
89 !> Get a list of species names
90 procedure :: get_spec_names
91 !> Get a species properties
92 procedure :: get_property_set
93 !> Get a species type
94 procedure :: get_type
95 !> Get a species phase
96 procedure :: get_phase
97 !> Get the absolute integration tolerance of a species
98 procedure :: get_abs_tol
99 !> Get a gas-phase species index in the \c
100 !! camp_camp_state::camp_state_t::state_var array. Note that
101 !! aerosol-phase species indices on the \c
102 !! camp_camp_state::camp_state_t::state_var array must be accessed from
103 !! \c camp_aero_rep_data::aero_rep_data_t::spec_state_id() for a
104 !! particular \ref camp_aero_rep "aerosol representation".
105 procedure :: gas_state_id
106 !> Get the name of a gas-phase species in the \c
107 !! camp_camp_state::camp_state_t::state_var array. Note that
108 !! aerosol-phase species names on the \c
109 !! camp_camp_state::camp_state_t::state_var array must be accessed from
110 !! \c camp_aero_rep_data::aero_rep_data_t::spec_state_id() for a
111 !! particular \ref camp_aero_rep "aerosol representation".
112 procedure :: gas_state_name
113 !> Print out the species data
114 procedure :: print => do_print
115 !> Finalize the chemical species data
116 final :: finalize
117
118 ! Private functions
119 !> Ensure there is enough room in the species dataset to add a
120 !! specified number of species
121 procedure, private :: ensure_size
122 !> Add a species
123 procedure, private :: add
124 !> Find a species index by name
125 procedure, private :: find
126 end type chem_spec_data_t
127
128 ! Constructor for chem_spec_data_t
129 interface chem_spec_data_t
130 procedure :: constructor
131 end interface chem_spec_data_t
132
133contains
134
135!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
136
137 !> Constructor for chem_spec_data_t
138 function constructor(init_size) result(new_obj)
139
140 !> A new set of chemical species
141 type(chem_spec_data_t), pointer :: new_obj
142 !> Number of species to allocate space for initially
143 integer(i_kind), intent(in), optional :: init_size
144
145 integer(i_kind) :: alloc_size
146
147 alloc_size = realloc_inc
148
149 if (present(init_size)) alloc_size = init_size
150 allocate(new_obj)
151 allocate(new_obj%spec_name(alloc_size))
152 allocate(new_obj%spec_type(alloc_size))
153 allocate(new_obj%spec_phase(alloc_size))
154 allocate(new_obj%property_set(alloc_size))
155
156 end function constructor
157
158!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
159
160 !> \page input_format_species Input JSON Object Format: Chemical Species
161 !!
162 !! A \c json object containing information about a \ref camp_species
163 !! "chemical species" has the following format:
164 !! \code{.json}
165 !! { "camp-data" : [
166 !! {
167 !! "name" : "my species name",
168 !! "type" : "CHEM_SPEC",
169 !! "phase" : "SPEC_PHASE",
170 !! "type" : "SPEC_TYPE",
171 !! "some property" : 123.34,
172 !! "some other property" : true,
173 !! "nested properties" : {
174 !! "sub prop 1" : 12.43,
175 !! "sub prop other" : "some text"
176 !! },
177 !! ...
178 !! },
179 !! {
180 !! "name" : "my species name",
181 !! "type" : "CHEM_SPEC",
182 !! "phase" : "SPEC_PHASE",
183 !! "type" : "SPEC_TYPE",
184 !! "some property" : 123.34,
185 !! ...
186 !! },
187 !! ...
188 !! ]}
189 !! \endcode
190 !! The key-value pair \b name is required and must contain the unique name
191 !! used for this species in the \ref input_format_mechanism
192 !! "mechanism object". (The same name cannot be used for a gas-phase species
193 !! and an aerosol-phase species.) The key-value pair \b type is also
194 !! required, and must be \b CHEM_SPEC.
195 !!
196 !! The key-value pair \b phase specifies the
197 !! phase in which the species exists and can be \b GAS or \b AEROSOL. When
198 !! the \b phase is not specified, it is assumed to be \b GAS. The \b type
199 !! can be \b VARIABLE, \b CONSTANT or \b PSSA. When a \b type is not
200 !! specified, it is assumed to be \b VARIABLE.
201 !!
202 !! All remaining data are optional and may include any valid \c json value,
203 !! including nested objects. Multilple entries with the same species name
204 !! will be merged into a single species, but duplicate property names for
205 !! the same species will cause an error. However, nested objects with the
206 !! same key name will be merged, if possible.
207
208!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
209
210 !> Load species from an input file
211#ifdef CAMP_USE_JSON
212 subroutine load(this, json, j_obj)
213
214 !> Species dataset
215 class(chem_spec_data_t), intent(inout) :: this
216 !> JSON core
217 type(json_core), pointer, intent(in) :: json
218 !> JSON object
219 type(json_value), pointer, intent(in) :: j_obj
220
221 type(json_value), pointer :: child, next
222 character(kind=json_ck, len=:), allocatable :: key, unicode_str_val
223 integer(kind=json_ik) :: var_type
224
225 character(len=:), allocatable :: spec_name, str_val
226 integer(kind=i_kind) :: spec_type, spec_phase
227 type(property_t), pointer :: property_set
228
229 ! allocate space for the species property set
230 property_set => property_t()
231
232 ! initialize type and phase
233 spec_type = chem_spec_unknown_type
234 spec_phase = chem_spec_unknown_phase
235
236 ! initialize species name in case of errors
237 spec_name = "unknown"
238
239 ! cycle through the species properties to find the name, type and phase
240 ! and load the remaining data into the species property set
241 next => null()
242 call json%get_child(j_obj, child)
243 do while (associated(child))
244 call json%info(child, name=key, var_type=var_type)
245
246 ! species name
247 if (key.eq."name") then
248 if (var_type.ne.json_string) call die_msg(181339359, &
249 "Received non-string species name")
250 call json%get(child, unicode_str_val)
251 spec_name = unicode_str_val
252
253 ! variable type
254 else if (key.eq."tracer type") then
255 if (var_type.ne.json_string) call die_msg(249541360, &
256 "Received non-string species type")
257 call json%get(child, unicode_str_val)
258 if (unicode_str_val.eq."VARIABLE") then
259 spec_type = chem_spec_variable
260 else if (unicode_str_val.eq."CONSTANT") then
261 spec_type = chem_spec_constant
262 else if (unicode_str_val.eq."PSSA") then
263 spec_type = chem_spec_pssa
264 else if (unicode_str_val.eq."ION_PAIR") then
265 spec_type = chem_spec_activity_coeff
266 spec_phase = chem_spec_aero_phase
267 else if (unicode_str_val.eq."ACTIVITY_COEFF") then
268 spec_type = chem_spec_activity_coeff
269 else
270 str_val = unicode_str_val
271 call die_msg(171550163, "Unknown chemical species type: "// &
272 str_val)
273 end if
274
275 ! species phase
276 else if (key.eq."phase") then
277 if (var_type.ne.json_string) call die_msg(273466657, &
278 "Received non-string species phase")
279 call json%get(child, unicode_str_val)
280 if (unicode_str_val.eq."GAS") then
281 call assert_msg(334314929, spec_type.ne.chem_spec_activity_coeff, &
282 "Activity coefficients and ion pairs cannot be gas-phase "// &
283 "species.")
284 spec_phase = chem_spec_gas_phase
285 else if (unicode_str_val.eq."AEROSOL") then
286 spec_phase = chem_spec_aero_phase
287 else
288 str_val = unicode_str_val
289 call die_msg(603704146, "Unknown chemical species phase: "// &
290 str_val)
291 end if
292
293 ! load remaining properties into the species property set
294 else if (key.ne."type") then
295 call property_set%load(json, child, .false., spec_name)
296 end if
297
298 call json%get_next(child, next)
299 child => next
300 end do
301
302 ! Add or update the species data
303 call this%add(spec_name, spec_type, spec_phase, property_set)
304
305 ! deallocate the temporary property set
306 deallocate(property_set)
307#else
308 subroutine load(this)
309
310 !> Species dataset
311 class(chem_spec_data_t), intent(inout) :: this
312
313 call warn_msg(627397948, "No support for input files")
314#endif
315 end subroutine load
316
317!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
318
319 !> Initialize the species set
320 subroutine initialize(this)
321
322 !> Species database
323 class(chem_spec_data_t), intent(inout) :: this
324
325 ! Species index
326 integer(kind=i_kind) :: i_spec
327
328 do i_spec = 1, this%num_spec
329
330 ! Set default value for type
331 if (this%spec_type(i_spec).eq.chem_spec_unknown_type) &
332 this%spec_type(i_spec) = chem_spec_variable
333
334 ! Set default value for phase
335 if (this%spec_phase(i_spec).eq.chem_spec_unknown_phase) &
336 this%spec_phase(i_spec) = chem_spec_gas_phase
337
338 end do
339
340 end subroutine initialize
341
342!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
343
344 !> Get the number of species with the given properties. If no properties
345 !! are specified, return the total number of species.
346 integer(kind=i_kind) function get_size(this, spec_type, spec_phase) &
347 result(num_spec)
348
349 !> Species database
350 class(chem_spec_data_t), intent(in) :: this
351 !> State variable type for the species
352 integer(kind=i_kind), intent(in), optional :: spec_type
353 !> Phase of the species
354 integer(kind=i_kind), intent(in), optional :: spec_phase
355
356 ! species index
357 integer(kind=i_kind) :: i_spec
358
359 ! add up the number of species
360 if (present(spec_type).or.present(spec_phase)) then
361 num_spec = 0
362 do i_spec = 1, this%num_spec
363 if (present(spec_type)) then
364 if (spec_type.ne.this%spec_type(i_spec)) cycle
365 end if
366 if (present(spec_phase)) then
367 if (spec_phase.ne.this%spec_phase(i_spec)) cycle
368 end if
369 num_spec = num_spec + 1
370 end do
371 else
372 num_spec = this%num_spec
373 end if
374
375 end function get_size
376
377!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
378
379 !> Check if a species name is in the set of chemical species
380 logical function exists(this, spec_name) result (found)
381
382 !> Species dataset
383 class(chem_spec_data_t), intent(in) :: this
384 !> Species name
385 character(len=*), intent(in) :: spec_name
386
387 ! Index of species
388 integer(kind=i_kind) :: i_spec
389
390 found = this%find(spec_name, i_spec)
391
392 end function exists
393
394!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
395
396 !> Get a list of species names
397 !!
398 !! If a type or phase are specified, only species with the specified values
399 !! are returned. Otherwise, all species names are returned.
400 function get_spec_names(this, spec_type, spec_phase) result (spec_names)
401
402 !> Species names
403 type(string_t), allocatable :: spec_names(:)
404 !> Species dataset
405 class(chem_spec_data_t), intent(in) :: this
406 !> State variable type for the species
407 integer(kind=i_kind), intent(in), optional :: spec_type
408 !> Phase of the species
409 integer(kind=i_kind), intent(in), optional :: spec_phase
410
411 ! species index and counter
412 integer(kind=i_kind) :: i_spec, num_spec
413
414 ! add up the number of species
415 if (present(spec_type).or.present(spec_phase)) then
416 num_spec = 0
417 do i_spec = 1, this%num_spec
418 if (present(spec_type)) then
419 if (spec_type.ne.this%spec_type(i_spec)) cycle
420 end if
421 if (present(spec_phase)) then
422 if (spec_phase.ne.this%spec_phase(i_spec)) cycle
423 end if
424 num_spec = num_spec + 1
425 end do
426 else
427 num_spec = this%num_spec
428 end if
429
430 ! Allocate space for the list
431 allocate(spec_names(num_spec))
432
433 ! Add species names to the list
434 num_spec = 0
435 do i_spec = 1, this%num_spec
436 if (present(spec_type)) then
437 if (spec_type.ne.this%spec_type(i_spec)) cycle
438 end if
439 if (present(spec_phase)) then
440 if (spec_phase.ne.this%spec_phase(i_spec)) cycle
441 end if
442 num_spec = num_spec + 1
443 spec_names(num_spec)%string = this%spec_name(i_spec)%string
444 end do
445
446 end function get_spec_names
447
448!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
449
450 !> Get a species property set. Returns true if the species is found, or
451 !! false otherwise.
452 logical function get_property_set(this, spec_name, property_set) &
453 result(found)
454
455 !> Species dataset
456 class(chem_spec_data_t), intent(in) :: this
457 !> Species name to find properties of
458 character(len=*), intent(in) :: spec_name
459 !> Pointer to species properties
460 type(property_t), pointer, intent(out) :: property_set
461
462 integer(i_kind) :: spec_id
463
464 property_set => null()
465 found = this%find(spec_name, spec_id)
466 if (found) property_set => this%property_set(spec_id)
467
468 end function get_property_set
469
470!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
471
472 !> Get a species type by species name. Returns true if the species is found
473 !! or false otherwise.
474 logical function get_type(this, spec_name, spec_type) result (found)
475
476 !> Species dataset
477 class(chem_spec_data_t), intent(in) :: this
478 !> Species name to find properties of
479 character(len=*), intent(in) :: spec_name
480 !> Species type
481 integer(kind=i_kind), intent(out) :: spec_type
482
483 integer(i_kind) :: spec_id
484
485 spec_type = chem_spec_unknown_type
486 found = this%find(spec_name, spec_id)
487 if (found) spec_type = this%spec_type(spec_id)
488
489 end function get_type
490
491!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
492
493 !> Get a species phase by name. Returns true if the species is found, or
494 !! false otherwise.
495 logical function get_phase(this, spec_name, spec_phase) &
496 result(found)
497
498 !> Species dataset
499 class(chem_spec_data_t), intent(in) :: this
500 !> Species name to find properties of
501 character(len=*), intent(in) :: spec_name
502 !> Species phase
503 integer(kind=i_kind), intent(out) :: spec_phase
504
505 integer(i_kind) :: spec_id
506
507 spec_phase = chem_spec_unknown_phase
508 found = this%find(spec_name, spec_id)
509 if (found) spec_phase = this%spec_phase(spec_id)
510
511 end function get_phase
512
513!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
514
515 !> Get the absolute integration tolerance of a species by name. Returns
516 !! true if the species is found, or false otherwise.
517 logical function get_abs_tol(this, spec_name, abs_tol) &
518 result(found)
519
520 !> Species dataset
521 class(chem_spec_data_t), intent(in) :: this
522 !> Species name
523 character(len=*), intent(in) :: spec_name
524 !> Absolute integration tolerance
525 real(kind=dp), intent(out) :: abs_tol
526
527 character(len=:), allocatable :: key
528 integer(kind=i_kind) :: spec_id
529 real(kind=dp) :: val
530
531 abs_tol = default_abs_tol
532 key = "absolute integration tolerance"
533 found = this%find(spec_name, spec_id)
534 if (found) then
535 if (this%property_set(spec_id)%get_real(key, val)) then
536 abs_tol = val
537 end if
538 end if
539
540 end function get_abs_tol
541
542!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
543
544 !> Get a gas-phase species index in the \c
545 !! camp_camp_state::camp_state_t::state_var array. Note that
546 !! aerosol-phase species indices on the \c
547 !! camp_camp_state::camp_state_t::state_var array must be accessed from
548 !! \c camp_aero_rep_data::aero_rep_data_t::spec_state_id() for a
549 !! particular \ref camp_aero_rep "aerosol representation". Returns a valid
550 !! state array index if the species is found, or 0 otherwise
551 integer(kind=i_kind) function gas_state_id(this, spec_name)
552
553 !> Species dataset
554 class(chem_spec_data_t), intent(in) :: this
555 !> Species name
556 character(len=*), intent(in) :: spec_name
557
558 integer(kind=i_kind) :: i_spec
559
560 gas_state_id = 0
561 do i_spec = 1, this%num_spec
562 if (this%spec_phase(i_spec).eq.chem_spec_gas_phase) then
564 if (trim(spec_name).eq.trim(this%spec_name(i_spec)%string)) return
565 end if
566 end do
567 gas_state_id = 0
568
569 end function gas_state_id
570
571!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
572
573 !> Get a gas-phase species name in the \c
574 !! camp_camp_state::camp_state_t::state_var array. Note that
575 !! aerosol-phase species names on the \c
576 !! camp_camp_state::camp_state_t::state_var array must be accessed from
577 !! \c camp_aero_rep_data::aero_rep_data_t::spec_state_id() for a
578 !! particular \ref camp_aero_rep "aerosol representation". Returns a valid
579 !! state array index if the species is found, or 0 otherwise
580 function gas_state_name(this, spec_id) result(spec_name)
581
582 !> Species name
583 character(len=:), allocatable :: spec_name
584 !> Species dataset
585 class(chem_spec_data_t), intent(in) :: this
586 !> Species id
587 integer(kind=i_kind), intent(in) :: spec_id
588
589 integer(kind=i_kind) :: gas_state_id, i_spec
590
591 gas_state_id = 0
592 do i_spec = 1, this%num_spec
593 if (this%spec_phase(i_spec).eq.chem_spec_gas_phase) &
595 if (gas_state_id.eq.spec_id) then
596 spec_name = trim(this%spec_name(i_spec)%string)
597 return
598 end if
599 end do
600 spec_name = ""
601
602 end function gas_state_name
603
604!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
605
606 !> Print out the species data
607 subroutine do_print(this, file_unit)
608
609 !> Chemical species data
610 class(chem_spec_data_t), intent(in) :: this
611 !> File unit for output
612 integer(kind=i_kind), optional, intent(in) :: file_unit
613
614 integer(kind=i_kind) :: i_spec
615 integer(kind=i_kind) :: f_unit
616 character(len=:), allocatable :: spec_phase, spec_type
617
618 f_unit = 6
619 if (present(file_unit)) f_unit = file_unit
620
621 write(f_unit,*) "Number of unique chemical species: ", this%num_spec
622 write(f_unit,*) "(The number of state variables may differ when a ", &
623 "species is present in more than one aerosol phase.)"
624 do i_spec = 1, this%num_spec
625 select case (this%spec_type(i_spec))
627 spec_type = "unknown"
628 case (chem_spec_variable)
629 spec_type = "variable"
630 case (chem_spec_constant)
631 spec_type = "constant"
632 case (chem_spec_pssa)
633 spec_type = "PSSA"
635 spec_type = "binary activity coefficient"
636 case default
637 spec_type = "invalid"
638 end select
639 select case (this%spec_phase(i_spec))
641 spec_phase = "unknown"
643 spec_phase = "gas"
645 spec_phase = "aerosol"
646 case default
647 spec_phase = "invalid"
648 end select
649
650 write(f_unit,*) " ", this%spec_name(i_spec)%string
651 write(f_unit,*) " phase: ", spec_phase, "; type: ", spec_type
652 call this%property_set(i_spec)%print(f_unit)
653 end do
654
655 end subroutine do_print
656
657!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
658
659 !> Finalize the chemical species data
660 elemental subroutine finalize(this)
661
662 !> Species dataset
663 type(chem_spec_data_t), intent(inout) :: this
664
665 deallocate(this%spec_name)
666 deallocate(this%spec_type)
667 deallocate(this%spec_phase)
668 deallocate(this%property_set)
669
670 end subroutine finalize
671
672!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
673
674 !> Ensure there is enough room in the species dataset to add a specified
675 !! number of species
676 subroutine ensure_size(this, num_spec)
677
678 !> Species dataset
679 class(chem_spec_data_t), intent(inout) :: this
680 !> Number of new species to ensure space for
681 integer(i_kind), intent(in) :: num_spec
682
683 integer :: new_size
684 type(string_t), pointer :: new_name(:)
685 integer(kind=i_kind), pointer :: new_type(:)
686 integer(kind=i_kind), pointer :: new_phase(:)
687 type(property_t), pointer :: new_property_set(:)
688
689 if (size(this%spec_name) .ge. this%num_spec + num_spec) return
690 new_size = this%num_spec + num_spec + realloc_inc
691 allocate(new_name(new_size))
692 allocate(new_type(new_size))
693 allocate(new_phase(new_size))
694 allocate(new_property_set(new_size))
695 new_name(1:this%num_spec) = this%spec_name(1:this%num_spec)
696 new_type(1:this%num_spec) = this%spec_type(1:this%num_spec)
697 new_phase(1:this%num_spec) = this%spec_phase(1:this%num_spec)
698 call this%property_set(1:this%num_spec)%move(new_property_set( &
699 1:this%num_spec))
700 deallocate(this%spec_name)
701 deallocate(this%spec_type)
702 deallocate(this%spec_phase)
703 deallocate(this%property_set)
704 this%spec_name => new_name
705 this%spec_type => new_type
706 this%spec_phase => new_phase
707 this%property_set => new_property_set
708
709 end subroutine ensure_size
710
711!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
712
713 !> Add a new chemical species
714 subroutine add(this, spec_name, spec_type, spec_phase, property_set)
715
716 !> Species dataset
717 class(chem_spec_data_t), intent(inout) :: this
718 !> Name of species to add
719 character(len=*), intent(in) :: spec_name
720 !> State variable type
721 integer(kind=i_kind), intent(inout) :: spec_type
722 !> Species phase
723 integer(kind=i_kind), intent(inout) :: spec_phase
724 !> Property set for new species
725 type(property_t), intent(inout), optional :: property_set
726
727 integer(kind=i_kind) :: i_spec
728
729 ! if the species exists, append the new data
730 if (this%find(spec_name, i_spec)) then
731
732 ! Check for a type mismatch
733 if (spec_type.eq.chem_spec_unknown_type) &
734 spec_type = this%spec_type(i_spec)
735 call assert_msg(596247182, &
736 this%spec_type(i_spec).eq.chem_spec_unknown_type &
737 .or.spec_type.eq.this%spec_type(i_spec), &
738 "Type mismatch for species "//spec_name)
739
740 ! Check for a phase mismatch
741 if (spec_phase.eq.chem_spec_unknown_phase) &
742 spec_phase = this%spec_phase(i_spec)
743 call assert_msg(612991075, &
744 this%spec_phase(i_spec).eq.chem_spec_unknown_phase &
745 .or.spec_phase.eq.this%spec_phase(i_spec), &
746 "Phase mismatch for species "//spec_name)
747
748 ! Update the species properties
749 this%spec_type(i_spec) = spec_type
750 this%spec_phase(i_spec) = spec_phase
751 call this%property_set(i_spec)%update(property_set, spec_name)
752
753 ! ... otherwise, create a new species
754 else
755 call this%ensure_size(1)
756 this%num_spec = this%num_spec + 1
757 this%spec_name(this%num_spec)%string = spec_name
758 this%spec_type(this%num_spec) = spec_type
759 this%spec_phase(this%num_spec) = spec_phase
760 if (present(property_set)) then
761 call property_set%move(this%property_set(this%num_spec))
762 else
763 this%property_set(this%num_spec) = property_t()
764 end if
765 end if
766
767 end subroutine add
768
769!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
770
771 !> Get the index of a chemical species by name. Returns true if the species
772 !! is found or false otherwise.
773 logical function find(this, spec_name, spec_id)
774
775 !> Species dataset
776 class(chem_spec_data_t), intent(in) :: this
777 !> Species name
778 character(len=*), intent(in) :: spec_name
779 !> Species id
780 integer(kind=i_kind), intent(out) :: spec_id
781
782 find = .true.
783 do spec_id = 1, this%num_spec
784 if (this%spec_name(spec_id)%string .eq. spec_name) then
785 return
786 end if
787 end do
788 find = .false.
789 spec_id = 0
790
791 end function find
792
793!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
794
795end module camp_chem_spec_data
The chem_spec_data_t structure and associated subroutines.
integer(kind=i_kind), parameter realloc_inc
Reallocation increment.
integer(kind=i_kind), parameter, public chem_spec_activity_coeff
logical function get_type(this, spec_name, spec_type)
Get a species type by species name. Returns true if the species is found or false otherwise.
logical function get_abs_tol(this, spec_name, abs_tol)
Get the absolute integration tolerance of a species by name. Returns true if the species is found,...
logical function find(this, spec_name, spec_id)
Get the index of a chemical species by name. Returns true if the species is found or false otherwise.
integer(kind=i_kind), parameter, public chem_spec_aero_phase
integer(kind=i_kind) function get_size(this, spec_type, spec_phase)
Get the number of species with the given properties. If no properties are specified,...
subroutine load(this, json, j_obj)
Load species from an input file.
integer(kind=i_kind), parameter, public chem_spec_unknown_type
State variable types (Must match values in camp_solver.c)
subroutine do_print(this, file_unit)
Print out the species data.
integer(kind=i_kind), parameter, public chem_spec_constant
logical function get_phase(this, spec_name, spec_phase)
Get a species phase by name. Returns true if the species is found, or false otherwise.
integer(kind=i_kind), parameter, public chem_spec_gas_phase
subroutine initialize(this)
Initialize the species set.
type(string_t) function, dimension(:), allocatable get_spec_names(this, spec_type, spec_phase)
Get a list of species names.
character(len=:) function, allocatable gas_state_name(this, spec_id)
Get a gas-phase species name in the camp_camp_state::camp_state_t::state_var array....
subroutine ensure_size(this, num_spec)
Ensure there is enough room in the species dataset to add a specified number of species.
type(chem_spec_data_t) function, pointer constructor(init_size)
Constructor for chem_spec_data_t.
integer(kind=i_kind), parameter, public chem_spec_unknown_phase
Species phase.
integer(kind=i_kind), parameter, public chem_spec_variable
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.
logical function exists(this, spec_name)
Check if a species name is in the set of chemical species.
integer(kind=i_kind) function gas_state_id(this, spec_name)
Get a gas-phase species index in the camp_camp_state::camp_state_t::state_var array....
elemental subroutine finalize(this)
Finalize the chemical species data.
integer(kind=i_kind), parameter, public chem_spec_pssa
subroutine add(this, spec_name, spec_type, spec_phase, property_set)
Add a new chemical species.
real(kind=dp), parameter default_abs_tol
Default absolute integration tolerance.
Physical constants.
Definition constants.F90:9
integer, parameter dp
Kind of a double precision real number.
Definition constants.F90:16
integer, parameter i_kind
Kind of an integer.
Definition constants.F90:21
The property_t structure and associated subroutines.
Definition property.F90:9
Common utility subroutines.
Definition util.F90:9
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_msg(code, warning_msg, already_warned)
Prints a warning message.
Definition util.F90:90
String type for building arrays of string of various size.
Definition util.F90:38