71 integer(kind=i_kind) :: num_spec = 0
73 type(
string_t),
pointer :: spec_name(:) => null()
75 integer(kind=i_kind),
pointer :: spec_type(:) => null()
77 integer(kind=i_kind),
pointer :: spec_phase(:) => null()
143 integer(i_kind),
intent(in),
optional :: init_size
145 integer(i_kind) :: alloc_size
149 if (
present(init_size)) alloc_size = init_size
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))
212 subroutine load(this, json, j_obj)
217 type(json_core),
pointer,
intent(in) :: json
219 type(json_value),
pointer,
intent(in) :: j_obj
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
225 character(len=:),
allocatable :: spec_name, str_val
226 integer(kind=i_kind) :: spec_type, spec_phase
237 spec_name =
"unknown"
242 call json%get_child(j_obj, child)
243 do while (
associated(child))
244 call json%info(child, name=key, var_type=var_type)
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
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
260 else if (unicode_str_val.eq.
"CONSTANT")
then
262 else if (unicode_str_val.eq.
"PSSA")
then
264 else if (unicode_str_val.eq.
"ION_PAIR")
then
267 else if (unicode_str_val.eq.
"ACTIVITY_COEFF")
then
270 str_val = unicode_str_val
271 call die_msg(171550163,
"Unknown chemical species type: "// &
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
282 "Activity coefficients and ion pairs cannot be gas-phase "// &
285 else if (unicode_str_val.eq.
"AEROSOL")
then
288 str_val = unicode_str_val
289 call die_msg(603704146,
"Unknown chemical species phase: "// &
294 else if (key.ne.
"type")
then
295 call property_set%load(json, child, .false., spec_name)
298 call json%get_next(child, next)
303 call this%add(spec_name, spec_type, spec_phase, property_set)
306 deallocate(property_set)
308 subroutine load(this)
313 call warn_msg(627397948,
"No support for input files")
326 integer(kind=i_kind) :: i_spec
328 do i_spec = 1, this%num_spec
346 integer(kind=i_kind) function get_size(this, spec_type, spec_phase) &
352 integer(kind=i_kind),
intent(in),
optional :: spec_type
354 integer(kind=i_kind),
intent(in),
optional :: spec_phase
357 integer(kind=i_kind) :: i_spec
360 if (
present(spec_type).or.
present(spec_phase))
then
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
366 if (
present(spec_phase))
then
367 if (spec_phase.ne.this%spec_phase(i_spec)) cycle
369 num_spec = num_spec + 1
372 num_spec = this%num_spec
380 logical function exists(this, spec_name)
result (found)
385 character(len=*),
intent(in) :: spec_name
388 integer(kind=i_kind) :: i_spec
390 found = this%find(spec_name, i_spec)
403 type(
string_t),
allocatable :: spec_names(:)
407 integer(kind=i_kind),
intent(in),
optional :: spec_type
409 integer(kind=i_kind),
intent(in),
optional :: spec_phase
412 integer(kind=i_kind) :: i_spec, num_spec
415 if (
present(spec_type).or.
present(spec_phase))
then
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
421 if (
present(spec_phase))
then
422 if (spec_phase.ne.this%spec_phase(i_spec)) cycle
424 num_spec = num_spec + 1
427 num_spec = this%num_spec
431 allocate(spec_names(num_spec))
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
439 if (
present(spec_phase))
then
440 if (spec_phase.ne.this%spec_phase(i_spec)) cycle
442 num_spec = num_spec + 1
443 spec_names(num_spec)%string = this%spec_name(i_spec)%string
458 character(len=*),
intent(in) :: spec_name
460 type(
property_t),
pointer,
intent(out) :: property_set
462 integer(i_kind) :: spec_id
464 property_set => null()
465 found = this%find(spec_name, spec_id)
466 if (found) property_set => this%property_set(spec_id)
474 logical function get_type(this, spec_name, spec_type)
result (found)
479 character(len=*),
intent(in) :: spec_name
481 integer(kind=i_kind),
intent(out) :: spec_type
483 integer(i_kind) :: spec_id
486 found = this%find(spec_name, spec_id)
487 if (found) spec_type = this%spec_type(spec_id)
495 logical function get_phase(this, spec_name, spec_phase) &
501 character(len=*),
intent(in) :: spec_name
503 integer(kind=i_kind),
intent(out) :: spec_phase
505 integer(i_kind) :: spec_id
508 found = this%find(spec_name, spec_id)
509 if (found) spec_phase = this%spec_phase(spec_id)
523 character(len=*),
intent(in) :: spec_name
525 real(kind=
dp),
intent(out) :: abs_tol
527 character(len=:),
allocatable :: key
528 integer(kind=i_kind) :: spec_id
532 key =
"absolute integration tolerance"
533 found = this%find(spec_name, spec_id)
535 if (this%property_set(spec_id)%get_real(key, val))
then
556 character(len=*),
intent(in) :: spec_name
558 integer(kind=i_kind) :: i_spec
561 do i_spec = 1, this%num_spec
564 if (trim(spec_name).eq.trim(this%spec_name(i_spec)%string))
return
583 character(len=:),
allocatable :: spec_name
587 integer(kind=i_kind),
intent(in) :: spec_id
592 do i_spec = 1, this%num_spec
596 spec_name = trim(this%spec_name(i_spec)%string)
612 integer(kind=i_kind),
optional,
intent(in) :: file_unit
614 integer(kind=i_kind) :: i_spec
615 integer(kind=i_kind) :: f_unit
616 character(len=:),
allocatable :: spec_phase, spec_type
619 if (
present(file_unit)) f_unit = file_unit
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"
629 spec_type =
"variable"
631 spec_type =
"constant"
635 spec_type =
"binary activity coefficient"
637 spec_type =
"invalid"
639 select case (this%spec_phase(i_spec))
641 spec_phase =
"unknown"
645 spec_phase =
"aerosol"
647 spec_phase =
"invalid"
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)
665 deallocate(this%spec_name)
666 deallocate(this%spec_type)
667 deallocate(this%spec_phase)
668 deallocate(this%property_set)
681 integer(i_kind),
intent(in) :: num_spec
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(:)
689 if (
size(this%spec_name) .ge. this%num_spec + num_spec)
return
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( &
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
714 subroutine add(this, spec_name, spec_type, spec_phase, property_set)
719 character(len=*),
intent(in) :: spec_name
721 integer(kind=i_kind),
intent(inout) :: spec_type
723 integer(kind=i_kind),
intent(inout) :: spec_phase
725 type(
property_t),
intent(inout),
optional :: property_set
727 integer(kind=i_kind) :: i_spec
730 if (this%find(spec_name, i_spec))
then
734 spec_type = this%spec_type(i_spec)
737 .or.spec_type.eq.this%spec_type(i_spec), &
738 "Type mismatch for species "//spec_name)
742 spec_phase = this%spec_phase(i_spec)
745 .or.spec_phase.eq.this%spec_phase(i_spec), &
746 "Phase mismatch for species "//spec_name)
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)
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))
763 this%property_set(this%num_spec) =
property_t()
773 logical function find(this, spec_name, spec_id)
778 character(len=*),
intent(in) :: spec_name
780 integer(kind=i_kind),
intent(out) :: spec_id
783 do spec_id = 1, this%num_spec
784 if (this%spec_name(spec_id)%string .eq. spec_name)
then
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.
integer, parameter dp
Kind of a double precision real number.
integer, parameter i_kind
Kind of an integer.
The property_t structure and associated subroutines.
Common utility subroutines.
subroutine die_msg(code, error_msg)
Error immediately.
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
subroutine warn_msg(code, warning_msg, already_warned)
Prints a warning message.
String type for building arrays of string of various size.