CAMP 1.0.0
Chemistry Across Multiple Phases
aero_phase_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_aero_phase_data module.
7
8!> \page camp_aero_phase CAMP: Aerosol Phase
9!!
10!! An \c camp_aero_phase_data::aero_phase_data_t object describes a distinct
11!! chemical phase within an aerosol. It is designed to allow the
12!! implementation of the chemical and mass transfer processes to be
13!! independent of the particular \ref camp_aero_rep "aerosol representation"
14!! used (e.g., bins, modes, single particles).
15!!
16!! A single \ref camp_aero_phase "aerosol phase" may be present in several
17!! \ref camp_aero_rep "aerosol representations" (e.g., an aqueous phase in a
18!! binned and a single-particle representation), but the \ref camp_species
19!! "chemical species" associated with a particular phase are constant
20!! throughout the model run. Once loaded, \ref camp_aero_phase
21!! "aerosol phases" are made available to any \ref input_format_aero_rep
22!! "aerosol representations" that want to implement them.
23!! \ref camp_aero_rep "Aerosol representations" are able to specify which
24!! phases they implement and how many instances of that phase are present in
25!! the \ref camp_aero_rep "aerosol representation". For example, a binned
26!! representation with 10 bins may implement 10 aqueous phases and 10 organic
27!! phases, whereas a single particle representation with a concentric shell
28!! structure of 3 layers may implement 3 of each phase (assuming the chemistry
29!! is solved for each particle individually).
30!!
31!! The set of \ref camp_aero_phase "aerosol phases" is made available to the
32!! \ref camp_mechanism "mechanism(s)" during model intialization. Reactions
33!! in the chemical mechanism are able to specify, by name, the phase
34!! in which they take place, and which species in that phase are involved.
35!! (How they decide this is up to the particular \ref camp_rxn
36!! "reaction type".) Any physical aerosol parameters, such as the surface
37!! area between phases, the particle radius, or the number concentration,
38!! required by a chemical reaction will be provided by the \ref camp_aero_rep
39!! "aerosol representation" at run time.
40!!
41!! The input format for an aerosol phase can be found \ref
42!! input_format_aero_phase "here".
43
44!> The abstract aero_phase_data_t structure and associated subroutines.
46
47#ifdef CAMP_USE_JSON
48 use json_module
49#endif
50#ifdef CAMP_USE_MPI
51 use mpi
52#endif
54 use camp_constants, only : i_kind, dp
55 use camp_mpi
58 use camp_util, only : die_msg, string_t
59
60 implicit none
61 private
62
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)
69
71
72 !> Reallocation increment
73 integer(kind=i_kind), parameter :: realloc_inc = 50
74
75 !> Aerosol phase data type
76 !!
77 !! \ref camp_aero_phase "Aerosol phase" information.
79 private
80 !> Name of the aerosol phase
81 character(len=:), allocatable :: phase_name
82 !> Number of species in the phase
83 integer(kind=i_kind) :: num_spec = 0
84 !> Species names. These are species that are present in the aerosol
85 !! phase. These species must exist in the \c
86 !! camp_camp_core::camp_core_t::chem_spec_data variable during
87 !! initialization.
88 type(string_t), pointer :: spec_name(:) => null()
89 !> Aerosol phase parameters. These will be available during
90 !! initialization, but not during solving.
91 type(property_t), pointer :: property_set => null()
92 !> Condensed phase data. Theses arrays will be available during
93 !! solving, and should contain any information required by the
94 !! functions of the aerosol phase that cannot be obtained
95 !! from the \c camp_camp_state::camp_state_t object. (floating-point)
96 real(kind=dp), allocatable, public :: condensed_data_real(:)
97 !> Condensed phase data. Theses arrays will be available during
98 !! solving, and should contain any information required by the
99 !! functions of the aerosol phase that cannot be obtained
100 !! from the \c camp_camp_state::camp_state_t object. (integer)
101 integer(kind=i_kind), allocatable, public :: condensed_data_int(:)
102 !> Pointer to the set of chemical species data
103 type(chem_spec_data_t), pointer :: chem_spec_data
104 contains
105 !> Load data from an input file
106 procedure :: load
107 !> Aerosol phase initialization
108 procedure :: initialize
109 !> Get the name of the aerosol phase
110 procedure :: name => get_name
111 !> Get the number of species in the phase
112 procedure :: size => get_size
113 !> Get the number of Jacobian row elements needed during solving
114 procedure :: num_jac_elem
115 !> Get property data associated with this phase
116 procedure :: get_property_set
117 !> Get a list of species names in this phase
118 procedure :: get_species_names
119 !> Get a species type by name
120 procedure :: get_species_type
121 !> Determine the number of bytes required to pack the given value
122 procedure :: pack_size
123 !> Packs the given value into the buffer, advancing position
124 procedure :: bin_pack
125 !> Unpacks the given value from the buffer, advancing position
126 procedure :: bin_unpack
127 !> Print the aerosol phase data
128 procedure :: print => do_print
129 !> Finalize the aerosol phase data
130 final :: finalize
131
132 ! Private functions
133 !> Ensure there is enough room in the species dataset to add a specified
134 !! number of species
135 procedure, private :: ensure_size
136 !> Add a species
137 procedure, private :: add
138 !> Find a species index by name
139 procedure, private :: find
140 end type aero_phase_data_t
141
142 ! Constructor for aero_phase_data_t
143 interface aero_phase_data_t
144 procedure :: constructor
145 end interface aero_phase_data_t
146
147 !> Pointer type for building arrays
149 type(aero_phase_data_t), pointer :: val => null()
150 contains
151 !> Dereference the pointer
152 procedure :: dereference
153 !> Finalize the pointer
155 end type aero_phase_data_ptr
156
157contains
158
159!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
160
161 !> Constructor for aero_phase_data_t
162 function constructor(phase_name, init_size) result(new_obj)
163
164 !> A new set of aerosol-phase species
165 type(aero_phase_data_t), pointer :: new_obj
166 !> Name of the aerosol phase
167 character(len=*), intent(in), optional :: phase_name
168 !> Number of species to allocate space for initially
169 integer(kind=i_kind), intent(in), optional :: init_size
170
171 integer(kind=i_kind) :: alloc_size = realloc_inc
172
173 if (present(init_size)) alloc_size = init_size
174 allocate(new_obj)
175 if (present(phase_name)) then
176 new_obj%phase_name = phase_name
177 else
178 new_obj%phase_name = ""
179 endif
180 allocate(new_obj%spec_name(alloc_size))
181
182 end function constructor
183
184!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
185
186 !> \page input_format_aero_phase Input JSON Object Format: Aerosol Phase
187 !!
188 !! A \c json object containing information about an \ref camp_aero_phase
189 !! "aerosol phase" has the following format:
190 !! \code{.json}
191 !! { "camp-data" : [
192 !! {
193 !! "name" : "my aerosol phase"
194 !! "type" : "AERO_PHASE"
195 !! "species" : [
196 !! "a species",
197 !! "another species",
198 !! ...
199 !! ],
200 !! ...
201 !! },
202 !! ...
203 !! ]}
204 !! \endcode
205 !! The key-value pair \b name is required and must contain the unique name
206 !! used for this \ref camp_aero_phase "aerosol phase" in the \ref
207 !! input_format_mechanism "mechanism". The key-value pair \b type is also
208 !! required and its value must be \b AERO_PHASE.
209 !!
210 !! A list of species names should be included in a key-value pair named
211 !! \b species whose value is an array of species names. These names must
212 !! correspond to \ref input_format_species "chemcical species" names.
213 !! \ref input_format_species "Chemical species" included in the \b species
214 !! array must have a \b phase of \b AEROSOL and must include key value pairs
215 !! \b molecular \b weight \b [\b kg \b mol-1]
216 !! and \b density \b [\b kg \b m-3].
217 !!
218 !! All other data is optional and may include any valid \c json value.
219 !! Multiple entries with the same aerosol phase \b name will be merged into
220 !! a single phase, but duplicate property names for the same phase will
221 !! cause an error.
222
223!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
224
225 !> Load species from an input file
226#ifdef CAMP_USE_JSON
227 subroutine load(this, json, j_obj)
228
229 !> Aerosol phase data
230 class(aero_phase_data_t), intent(inout) :: this
231 !> JSON core
232 type(json_core), pointer, intent(in) :: json
233 !> JSON object
234 type(json_value), pointer, intent(in) :: j_obj
235
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
239
240 character(len=:), allocatable :: str_val
241 type(property_t), pointer :: property_set
242
243 ! allocate space for the phase property set
244 property_set => property_t()
245
246 ! cycle through the phase properties to find the name and species
247 ! and load the remaining data into the phase property set
248 next => null()
249 call json%get_child(j_obj, child)
250 do while (associated(child))
251 call json%info(child, name=key, var_type=var_type)
252
253 ! phase name
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
259
260 ! chemical species in the phase
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: "//&
264 to_string(var_type))
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)
274 species => next
275 end do
276
277 ! load remaining properties into the phase property set
278 else if (key.ne."type") then
279 call property_set%load(json, child, .false., this%phase_name)
280 end if
281
282 call json%get_next(child, next)
283 child => next
284 end do
285
286 ! save the property set
287 if (associated(this%property_set)) then
288 call this%property_set%update(property_set, this%phase_name)
289 deallocate (property_set)
290 else
291 this%property_set => property_set
292 end if
293#else
294 subroutine load(this)
295
296 !> Aerosol phase data
297 class(aero_phase_data_t), intent(in) :: this
298
299 call warn_msg(236665532, "No support for input files.")
300#endif
301 end subroutine load
302
303!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
304
305 !> Initialize the aerosol phase data, validating species names.
306 subroutine initialize(this, chem_spec_data)
307
308 !> Aerosol phase data
309 class(aero_phase_data_t), intent(inout) :: this
310 !> Chemical species data
311 type(chem_spec_data_t), target, intent(in) :: chem_spec_data
312
313 type(property_t), pointer :: spec_props
314 integer(kind=i_kind) :: i_spec, i_spec_phase_type
315 character(len=:), allocatable :: key_name
316
317 ! Allocate space for the condensed data arrays
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))
320
321 ! Set the number of species
322 num_state_var_ = this%num_spec
323
324 ! Find the aerosol-phase species, and save their needed properties
325 do i_spec = 1, num_state_var_
326
327 ! Get the species properties
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//"'")
333
334 ! Get the species type
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//"'")
340
341 ! Make sure the species is an aerosol-phase species
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)
346 call assert_msg(861388228, i_spec_phase_type.eq.chem_spec_aero_phase, &
347 "Trying to add non-aerosol phase species to aerosol phase "// &
348 this%phase_name//"; species: "//this%spec_name(i_spec)%string)
349
350 ! get the molecular weight and density of species
351 ! present in the phase
352 if (spec_type_(i_spec).eq.chem_spec_variable .or. &
353 spec_type_(i_spec).eq.chem_spec_constant .or. &
354 spec_type_(i_spec).eq.chem_spec_pssa) then
355
356 ! Get the molecular weight
357 key_name = "molecular weight [kg mol-1]"
358 call assert_msg(512254139, &
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//"'")
363
364 ! Get the density
365 key_name = "density [kg m-3]"
366 call assert_msg(224966878, &
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//"'")
371
372 ! activity coefficients do not need molecular weight or density
373 else
374
375 mw_(i_spec) = 0.0
376 density_(i_spec) = 0.0
377
378 end if
379
380 end do
381
382 ! Save a pointer to the chemical species data
383 this%chem_spec_data => chem_spec_data
384
385 end subroutine initialize
386
387!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
388
389 !> Get the aerosol phase name
390 function get_name(this) result (phase_name)
391
392 !> The name of the aerosol phase
393 character(len=:), allocatable :: phase_name
394 !> Aerosol phase data
395 class(aero_phase_data_t), intent(in) :: this
396
397 phase_name = this%phase_name
398
399 end function get_name
400
401!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
402
403 !> Get the number of species in the phase
404 integer(kind=i_kind) function get_size(this) result(num_spec)
405
406 !> Aerosol phase data
407 class(aero_phase_data_t), intent(in) :: this
408
409 num_spec = this%num_spec
410
411 end function get_size
412
413!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
414
415 !> Get the number of Jacobian row elements needed during solving
416 integer(kind=i_kind) function num_jac_elem(this) result(num_elem)
417
418 !> Aerosol phase data
419 class(aero_phase_data_t), intent(in) :: this
420
421 integer :: i_spec
422
423 num_elem = 0
424 do i_spec = 1, this%num_spec
425 if (spec_type_(i_spec) == chem_spec_variable .or. &
426 spec_type_(i_spec) == chem_spec_constant .or. &
427 spec_type_(i_spec) == chem_spec_pssa) then
428 num_elem = num_elem + 1
429 end if
430 end do
431
432 end function num_jac_elem
433
434!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
435
436 !> Get the aerosol phase property set
437 function get_property_set(this) result (property_set)
438
439 !> A pointer to the aerosol phase property set
440 class(property_t), pointer :: property_set
441 !> Aerosol phase data
442 class(aero_phase_data_t), intent(in) :: this
443
444 property_set => this%property_set
445
446 end function get_property_set
447
448!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
449
450 !> Get an aerosol phase species name
451 function get_species_names(this) result (spec_names)
452
453 !> Names of species in this phase
454 type(string_t), allocatable :: spec_names(:)
455 !> Aerosol phase data
456 class(aero_phase_data_t), intent(in) :: this
457
458 integer(kind=i_kind) :: i_spec
459
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
463 end do
464
465 end function get_species_names
466
467!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
468
469 !> Get an aerosol phase species type
470 function get_species_type(this, spec_name) result (spec_type)
471
472 !> The type of a species in this phase
473 integer(kind=i_kind) :: spec_type
474 !> Aerosol phase data
475 class(aero_phase_data_t), intent(in) :: this
476 !> Name of the species
477 character(len=*), intent(in) :: spec_name
478
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))
483
484 end function get_species_type
485
486!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
487
488 !> Determine the size of a binary required to pack the aerosol
489 !! representation data
490 integer(kind=i_kind) function pack_size(this, comm)
491
492 !> Aerosol representation data
493 class(aero_phase_data_t), intent(in) :: this
494 !> MPI communicator
495 integer, intent(in) :: comm
496
497 pack_size = &
498 camp_mpi_pack_size_real_array(this%condensed_data_real, comm) + &
499 camp_mpi_pack_size_integer_array(this%condensed_data_int, comm)
500
501 end function pack_size
502
503!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
504
505 !> Pack the given value to the buffer, advancing position
506 subroutine bin_pack(this, buffer, pos, comm)
507
508 !> Aerosol representation data
509 class(aero_phase_data_t), intent(in) :: this
510 !> Memory buffer
511 character, intent(inout) :: buffer(:)
512 !> Current buffer position
513 integer, intent(inout) :: pos
514 !> MPI communicator
515 integer, intent(in) :: comm
516
517#ifdef CAMP_USE_MPI
518 integer :: prev_position
519
520 prev_position = pos
521 call camp_mpi_pack_real_array(buffer, pos, this%condensed_data_real, comm)
522 call camp_mpi_pack_integer_array(buffer, pos, this%condensed_data_int,comm)
523 call assert(561436372, &
524 pos - prev_position <= this%pack_size(comm))
525#endif
526
527 end subroutine bin_pack
528
529!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
530
531 !> Unpack the given value from the buffer, advancing position
532 subroutine bin_unpack(this, buffer, pos, comm)
533
534 !> Aerosol representation data
535 class(aero_phase_data_t), intent(out) :: this
536 !> Memory buffer
537 character, intent(inout) :: buffer(:)
538 !> Current buffer position
539 integer, intent(inout) :: pos
540 !> MPI communicator
541 integer, intent(in) :: comm
542
543#ifdef CAMP_USE_MPI
544 integer :: prev_position
545
546 prev_position = pos
547 call camp_mpi_unpack_real_array(buffer, pos, this%condensed_data_real,comm)
548 call camp_mpi_unpack_integer_array(buffer, pos, this%condensed_data_int, &
549 comm)
550 call assert(219217030, &
551 pos - prev_position <= this%pack_size(comm))
552#endif
553
554 end subroutine bin_unpack
555
556!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
557
558 !> Print out the aerosol phase data
559 subroutine do_print(this, file_unit)
560
561 !> Aerosol phase data
562 class(aero_phase_data_t), intent(in) :: this
563 !> File unit for output
564 integer(kind=i_kind), optional :: file_unit
565
566 integer(kind=i_kind) :: f_unit
567 integer(kind=i_kind) :: i_spec
568
569 f_unit = 6
570
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
577 end do
578 write(f_unit,*) "]"
579 call this%property_set%print(f_unit)
580 write(f_unit,*) "End aerosol phase: ", this%phase_name
581
582 end subroutine do_print
583
584!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
585
586 !> Finalize the aerosol phase data
587 elemental subroutine finalize(this)
588
589 !> Aerosol phase data
590 type(aero_phase_data_t), intent(inout) :: this
591
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)
599
600 end subroutine finalize
601
602!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
603
604 !> Ensure there is enough room in the species dataset to add a specified
605 !! number of species
606 subroutine ensure_size(this, num_spec)
607
608 !> Aerosol phase data
609 class(aero_phase_data_t), intent(inout) :: this
610 !> Number of new species to ensure space for
611 integer(kind=i_kind), intent(in) :: num_spec
612
613 integer :: new_size
614 type(string_t), pointer :: new_name(:)
615
616 if (size(this%spec_name) .ge. this%num_spec + num_spec) return
617 new_size = this%num_spec + num_spec + realloc_inc
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
622
623 end subroutine ensure_size
624
625!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
626
627 !> Add a new chemical species to the phase
628 subroutine add(this, spec_name)
629
630 !> Aerosol phase data
631 class(aero_phase_data_t), intent(inout) :: this
632 !> Name of the species to add
633 character(len=*), intent(in) :: spec_name
634
635 integer(kind=i_kind) :: i_spec
636
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())
641 return
642 end if
643 call this%ensure_size(1)
644 this%num_spec = this%num_spec + 1
645 this%spec_name(this%num_spec)%string = spec_name
646
647 end subroutine add
648
649!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
650
651 !> Get the index of an aerosol-phase species by name. Return 0 if the
652 !! species is not found
653 integer(kind=i_kind) function find(this, spec_name) &
654 result(spec_id)
655
656 !> Aerosol phase data
657 class(aero_phase_data_t), intent(in) :: this
658 !> Species name
659 character(len=*), intent(in) :: spec_name
660
661 integer(kind=i_kind) :: i_spec
662
663 spec_id = 0
664 do i_spec = 1, this%num_spec
665 if (this%spec_name(i_spec)%string .eq. spec_name) then
666 spec_id = i_spec
667 return
668 end if
669 end do
670
671 end function find
672
673!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
674
675 !> Dereference a pointer to aerosol phase data
676 elemental subroutine dereference(this)
677
678 !> Pointer to aerosol phase data
679 class(aero_phase_data_ptr), intent(inout) :: this
680
681 this%val => null()
682
683 end subroutine dereference
684
685!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
686
687 !> Finalize a pointer to aerosol phase data
688 elemental subroutine ptr_finalize(this)
689
690 !> Pointer to aerosol phase data
691 type(aero_phase_data_ptr), intent(inout) :: this
692
693 if (associated(this%val)) deallocate(this%val)
694
695 end subroutine ptr_finalize
696
697!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
698
699#undef NUM_STATE_VAR_
700#undef NUM_INT_PROP_
701#undef NUM_REAL_PROP_
702#undef SPEC_TYPE_
703#undef MW_
704#undef DENSITY_
705
706end module camp_aero_phase_data
Interface for to_string functions.
Definition util.F90:32
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.
Definition camp_state.F90:9
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
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
Wrapper functions for MPI.
Definition mpi.F90:13
subroutine camp_mpi_pack_integer_array(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition mpi.F90:858
subroutine camp_mpi_pack_real_array(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition mpi.F90:899
subroutine camp_mpi_unpack_integer_array(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition mpi.F90:1201
integer function camp_mpi_pack_size_real_array(val, comm)
Determines the number of bytes required to pack the given value.
Definition mpi.F90:578
integer function camp_mpi_pack_size_integer_array(val, comm)
Determines the number of bytes required to pack the given value.
Definition mpi.F90:540
subroutine camp_mpi_unpack_real_array(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition mpi.F90:1242
The property_t structure and associated subroutines.
Definition property.F90:9
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_msg(code, warning_msg, already_warned)
Prints a warning message.
Definition util.F90:90
Pointer type for building arrays.
String type for building arrays of string of various size.
Definition util.F90:38