63#define NUM_STATE_VAR_ this%condensed_data_int(1)
64#define NUM_INT_PROP_ 1
65#define NUM_REAL_PROP_ 0
66#define SPEC_TYPE_(x) this%condensed_data_int(NUM_INT_PROP_+x)
67#define MW_(x) this%condensed_data_real(NUM_REAL_PROP_+x)
68#define DENSITY_(x) this%condensed_data_real(NUM_REAL_PROP_+NUM_STATE_VAR_+x)
81 character(len=:),
allocatable :: phase_name
83 integer(kind=i_kind) :: num_spec = 0
88 type(
string_t),
pointer :: spec_name(:) => null()
96 real(kind=
dp),
allocatable,
public :: condensed_data_real(:)
101 integer(kind=i_kind),
allocatable,
public :: condensed_data_int(:)
167 character(len=*),
intent(in),
optional :: phase_name
169 integer(kind=i_kind),
intent(in),
optional :: init_size
173 if (
present(init_size)) alloc_size = init_size
175 if (
present(phase_name))
then
176 new_obj%phase_name = phase_name
178 new_obj%phase_name =
""
180 allocate(new_obj%spec_name(alloc_size))
227 subroutine load(this, json, j_obj)
232 type(json_core),
pointer,
intent(in) :: json
234 type(json_value),
pointer,
intent(in) :: j_obj
236 type(json_value),
pointer :: child, next, species
237 character(kind=json_ck, len=:),
allocatable :: key, unicode_str_val
238 integer(kind=i_kind) :: var_type
240 character(len=:),
allocatable :: str_val
249 call json%get_child(j_obj, child)
250 do while (
associated(child))
251 call json%info(child, name=key, var_type=var_type)
254 if (key.eq.
"name")
then
255 if (var_type.ne.json_string)
call die_msg(429142134, &
256 "Received non-string aerosol phase name.")
257 call json%get(child, unicode_str_val)
258 this%phase_name = unicode_str_val
261 else if (key.eq.
"species")
then
262 if (var_type.ne.json_array)
call die_msg(293312378, &
263 "Received non-array list of aerosol phase species: "//&
265 call json%get_child(child, species)
266 do while (
associated(species))
267 call json%info(species, var_type=var_type)
268 if (var_type.ne.json_string)
call die_msg(669858868, &
269 "Received non-string aerosol phase species name.")
270 call json%get(species, unicode_str_val)
271 str_val = unicode_str_val
272 call this%add(str_val)
273 call json%get_next(species, next)
278 else if (key.ne.
"type")
then
279 call property_set%load(json, child, .false., this%phase_name)
282 call json%get_next(child, next)
287 if (
associated(this%property_set))
then
288 call this%property_set%update(property_set, this%phase_name)
289 deallocate (property_set)
291 this%property_set => property_set
294 subroutine load(this)
299 call warn_msg(236665532,
"No support for input files.")
314 integer(kind=i_kind) :: i_spec, i_spec_phase_type
315 character(len=:),
allocatable :: key_name
318 allocate(this%condensed_data_int(num_int_prop_+this%num_spec))
319 allocate(this%condensed_data_real(num_real_prop_+2*this%num_spec))
322 num_state_var_ = this%num_spec
325 do i_spec = 1, num_state_var_
328 call assert_msg(140971956, chem_spec_data%get_property_set( &
329 this%spec_name(i_spec)%string, spec_props), &
330 "Missing property set for species '"// &
331 this%spec_name(i_spec)%string// &
332 "' in aerosol phase '"//this%phase_name//
"'")
335 call assert_msg(129442398, chem_spec_data%get_type( &
336 this%spec_name(i_spec)%string, spec_type_(i_spec)), &
337 "Missing type for species '"// &
338 this%spec_name(i_spec)%string// &
339 "' in aerosol phase '"//this%phase_name//
"'")
342 call assert_msg(136357145, chem_spec_data%get_phase( &
343 this%spec_name(i_spec)%string, i_spec_phase_type ), &
344 "Error getting phase for species "// &
345 this%spec_name(i_spec)%string)
347 "Trying to add non-aerosol phase species to aerosol phase "// &
348 this%phase_name//
"; species: "//this%spec_name(i_spec)%string)
357 key_name =
"molecular weight [kg mol-1]"
359 spec_props%get_real(key_name, mw_(i_spec)), &
360 "Missing molecular weight for species '"// &
361 this%spec_name(i_spec)%string// &
362 "' in aerosol phase '"//this%phase_name//
"'")
365 key_name =
"density [kg m-3]"
367 spec_props%get_real(key_name, density_(i_spec)), &
368 "Missing density for species '"// &
369 this%spec_name(i_spec)%string// &
370 "' in aerosol phase '"//this%phase_name//
"'")
376 density_(i_spec) = 0.0
383 this%chem_spec_data => chem_spec_data
393 character(len=:),
allocatable :: phase_name
397 phase_name = this%phase_name
404 integer(kind=i_kind) function get_size(this)
result(num_spec)
409 num_spec = this%num_spec
424 do i_spec = 1, this%num_spec
428 num_elem = num_elem + 1
444 property_set => this%property_set
454 type(
string_t),
allocatable :: spec_names(:)
458 integer(kind=i_kind) :: i_spec
460 allocate(spec_names(this%num_spec))
461 do i_spec = 1, this%num_spec
462 spec_names(i_spec)%string = this%spec_name(i_spec)%string
473 integer(kind=i_kind) :: spec_type
477 character(len=*),
intent(in) :: spec_name
479 call assert_msg(163269315, this%find(spec_name).gt.0, &
480 "Species '"//spec_name//
"' is not in aerosol phase '"// &
481 this%phase_name//
"'.")
482 call assert(255786656, this%chem_spec_data%get_type(spec_name, spec_type))
495 integer,
intent(in) :: comm
511 character,
intent(inout) :: buffer(:)
513 integer,
intent(inout) :: pos
515 integer,
intent(in) :: comm
518 integer :: prev_position
524 pos - prev_position <= this%pack_size(comm))
537 character,
intent(inout) :: buffer(:)
539 integer,
intent(inout) :: pos
541 integer,
intent(in) :: comm
544 integer :: prev_position
551 pos - prev_position <= this%pack_size(comm))
564 integer(kind=i_kind),
optional :: file_unit
566 integer(kind=i_kind) :: f_unit
567 integer(kind=i_kind) :: i_spec
571 if (
present(file_unit)) f_unit = file_unit
572 write(f_unit,*)
"Aerosol phase: ", this%phase_name
573 write(f_unit,*)
"Number of species: ", this%num_spec
574 write(f_unit,*)
"Species: ["
575 do i_spec = 1, this%num_spec
576 write(*,*) this%spec_name(i_spec)%string
579 call this%property_set%print(f_unit)
580 write(f_unit,*)
"End aerosol phase: ", this%phase_name
592 if (
allocated(this%phase_name))
deallocate(this%phase_name)
593 if (
associated(this%spec_name))
deallocate(this%spec_name)
594 if (
associated(this%property_set))
deallocate(this%property_set)
595 if (
allocated(this%condensed_data_real)) &
596 deallocate(this%condensed_data_real)
597 if (
allocated(this%condensed_data_int)) &
598 deallocate(this%condensed_data_int)
611 integer(kind=i_kind),
intent(in) :: num_spec
614 type(
string_t),
pointer :: new_name(:)
616 if (
size(this%spec_name) .ge. this%num_spec + num_spec)
return
618 allocate(new_name(new_size))
619 new_name(1:this%num_spec) = this%spec_name(1:this%num_spec)
620 deallocate(this%spec_name)
621 this%spec_name => new_name
628 subroutine add(this, spec_name)
633 character(len=*),
intent(in) :: spec_name
635 integer(kind=i_kind) :: i_spec
637 i_spec = this%find(spec_name)
638 if (i_spec.ne.0)
then
639 call warn_msg(980242449,
"Species "//spec_name//&
640 " added more than once to phase "//this%name())
643 call this%ensure_size(1)
644 this%num_spec = this%num_spec + 1
645 this%spec_name(this%num_spec)%string = spec_name
653 integer(kind=i_kind) function find(this, spec_name) &
659 character(len=*),
intent(in) :: spec_name
661 integer(kind=i_kind) :: i_spec
664 do i_spec = 1, this%num_spec
665 if (this%spec_name(i_spec)%string .eq. spec_name)
then
693 if (
associated(this%val))
deallocate(this%val)
Interface for to_string functions.
The abstract aero_phase_data_t structure and associated subroutines.
integer(kind=i_kind), parameter realloc_inc
Reallocation increment.
elemental subroutine finalize(this)
Finalize the aerosol phase data.
subroutine bin_unpack(this, buffer, pos, comm)
Unpack the given value from the buffer, advancing position.
integer(kind=i_kind) function get_size(this)
Get the number of species in the phase.
subroutine initialize(this, chem_spec_data)
Initialize the aerosol phase data, validating species names.
class(property_t) function, pointer get_property_set(this)
Get the aerosol phase property set.
character(len=:) function, allocatable get_name(this)
Get the aerosol phase name.
type(string_t) function, dimension(:), allocatable get_species_names(this)
Get an aerosol phase species name.
type(aero_phase_data_t) function, pointer constructor(phase_name, init_size)
Constructor for aero_phase_data_t.
subroutine bin_pack(this, buffer, pos, comm)
Pack the given value to the buffer, advancing position.
subroutine load(this, json, j_obj)
Load species from an input file.
subroutine add(this, spec_name)
Add a new chemical species to the phase.
integer(kind=i_kind) function num_jac_elem(this)
Get the number of Jacobian row elements needed during solving.
subroutine ensure_size(this, num_spec)
Ensure there is enough room in the species dataset to add a specified number of species.
integer(kind=i_kind) function find(this, spec_name)
Get the index of an aerosol-phase species by name. Return 0 if the species is not found.
integer(kind=i_kind) function pack_size(this, comm)
Determine the size of a binary required to pack the aerosol representation data.
elemental subroutine ptr_finalize(this)
Finalize a pointer to aerosol phase data.
subroutine do_print(this, file_unit)
Print out the aerosol phase data.
integer(kind=i_kind) function get_species_type(this, spec_name)
Get an aerosol phase species type.
elemental subroutine dereference(this)
Dereference a pointer to aerosol phase data.
The camp_state_t structure and associated subroutines.
The chem_spec_data_t structure and associated subroutines.
integer(kind=i_kind), parameter, public chem_spec_aero_phase
integer(kind=i_kind), parameter, public chem_spec_constant
integer(kind=i_kind), parameter, public chem_spec_variable
integer(kind=i_kind), parameter, public chem_spec_pssa
integer, parameter dp
Kind of a double precision real number.
integer, parameter i_kind
Kind of an integer.
Wrapper functions for MPI.
subroutine camp_mpi_pack_integer_array(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
subroutine camp_mpi_pack_real_array(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
subroutine camp_mpi_unpack_integer_array(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
integer function camp_mpi_pack_size_real_array(val, comm)
Determines the number of bytes required to pack the given value.
integer function camp_mpi_pack_size_integer_array(val, comm)
Determines the number of bytes required to pack the given value.
subroutine camp_mpi_unpack_real_array(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
The property_t structure and associated subroutines.
Common utility subroutines.
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
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.
Pointer type for building arrays.
String type for building arrays of string of various size.