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.
78 type :: aero_phase_data_t
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 !!> Property set associated with species. When species
93 !! are input in JSON files as objects, associated properties
94 !! can be included.
95 type(property_ptr), allocatable :: spec_property_set(:)
96 !> Condensed phase data. Theses arrays will be available during
97 !! solving, and should contain any information required by the
98 !! functions of the aerosol phase that cannot be obtained
99 !! from the \c camp_camp_state::camp_state_t object. (floating-point)
100 real(kind=dp), allocatable, public :: condensed_data_real(:)
101 !> Condensed phase data. Theses arrays will be available during
102 !! solving, and should contain any information required by the
103 !! functions of the aerosol phase that cannot be obtained
104 !! from the \c camp_camp_state::camp_state_t object. (integer)
105 integer(kind=i_kind), allocatable, public :: condensed_data_int(:)
106 !> Pointer to the set of chemical species data
107 type(chem_spec_data_t), pointer :: chem_spec_data
108 contains
109 !> Load data from an input file
110 procedure :: load
111 !> Aerosol phase initialization
112 procedure :: initialize
113 !> Get the name of the aerosol phase
114 procedure :: name => get_name
115 !> Get the number of species in the phase
116 procedure :: size => get_size
117 !> Get the number of Jacobian row elements needed during solving
118 procedure :: num_jac_elem
119 !> Get property data associated with this phase
120 procedure :: get_property_set
121 !> Get property data associated with a species in a phase
123 !> Get a list of species names in this phase
124 procedure :: get_species_names
125 !> Get a species type by name
126 procedure :: get_species_type
127 !> Determine the number of bytes required to pack the given value
128 procedure :: pack_size
129 !> Packs the given value into the buffer, advancing position
130 procedure :: bin_pack
131 !> Unpacks the given value from the buffer, advancing position
132 procedure :: bin_unpack
133 !> Print the aerosol phase data
134 procedure :: print => do_print
135 !> Finalize the aerosol phase data
137
138 ! Private functions
139 !> Ensure there is enough room in the species dataset to add a specified
140 !! number of species
141 procedure, private :: ensure_size
142 !> Add a species
143 procedure, private :: add
144 !> Find a species index by name
145 procedure, private :: find
146 end type aero_phase_data_t
147
148 ! Constructor for aero_phase_data_t
150 procedure :: constructor
151 end interface aero_phase_data_t
152
153 !> Pointer type for building arrays
155 type(aero_phase_data_t), pointer :: val => null()
156 contains
157 !> Dereference the pointer
158 procedure :: dereference
159 !> Finalize the pointer
161 end type aero_phase_data_ptr
162
163contains
164
165!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
166
167 !> Constructor for aero_phase_data_t
168 function constructor(phase_name, init_size) result(new_obj)
169
170 !> A new set of aerosol-phase species
171 type(aero_phase_data_t), pointer :: new_obj
172 !> Name of the aerosol phase
173 character(len=*), intent(in), optional :: phase_name
174 !> Number of species to allocate space for initially
175 integer(kind=i_kind), intent(in), optional :: init_size
176
177 integer(kind=i_kind) :: alloc_size = realloc_inc
178
179 if (present(init_size)) alloc_size = init_size
180 allocate(new_obj)
181 if (present(phase_name)) then
182 new_obj%phase_name = phase_name
183 else
184 new_obj%phase_name = ""
185 endif
186 allocate(new_obj%spec_name(alloc_size))
187
188 end function constructor
189
190!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
191
192 !> \page input_format_aero_phase Input JSON Object Format: Aerosol Phase
193 !!
194 !! A \c json object containing information about an \ref camp_aero_phase
195 !! "aerosol phase" has the following format:
196 !! \code{.json}
197 !! { "camp-data" : [
198 !! {
199 !! "name" : "my aerosol phase"
200 !! "type" : "AERO_PHASE"
201 !! "species" : [
202 !! "a species",
203 !! "another species",
204 !! ...
205 !! ],
206 !! ...
207 !! },
208 !! ...
209 !! ]}
210 !! \endcode
211 !! The key-value pair \b name is required and must contain the unique name
212 !! used for this \ref camp_aero_phase "aerosol phase" in the \ref
213 !! input_format_mechanism "mechanism". The key-value pair \b type is also
214 !! required and its value must be \b AERO_PHASE.
215 !!
216 !! A list of species names should be included in a key-value pair named
217 !! \b species whose value is an array of species names. These names must
218 !! correspond to \ref input_format_species "chemcical species" names.
219 !! \ref input_format_species "Chemical species" included in the \b species
220 !! array must have a \b phase of \b AEROSOL and must include key value pairs
221 !! \b molecular \b weight \b [\b kg \b mol-1]
222 !! and \b density \b [\b kg \b m-3].
223 !!
224 !! All other data is optional and may include any valid \c json value.
225 !! Multiple entries with the same aerosol phase \b name will be merged into
226 !! a single phase, but duplicate property names for the same phase will
227 !! cause an error.
228
229!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
230
231 !> Load species from an input file
232#ifdef CAMP_USE_JSON
233 subroutine load(this, json, j_obj)
234
235 !> Aerosol phase data
236 class(aero_phase_data_t), intent(inout) :: this
237 !> JSON core
238 type(json_core), pointer, intent(in) :: json
239 !> JSON object
240 type(json_value), pointer, intent(in) :: j_obj
241
242 type(json_value), pointer :: child, next, species, species_child
243 type(json_value), pointer :: species_obj
244 character(kind=json_ck, len=:), allocatable :: key, unicode_str_val
245 character(kind=json_ck, len=:), allocatable :: species_name
246 integer(kind=i_kind) :: var_type, num_spec, i_spec
247
248 character(len=:), allocatable :: str_val
249 type(property_t), pointer :: property_set
250 type(property_t), pointer :: spec_property_set
251
252 ! allocate space for the phase property set
253 property_set => property_t()
254
255 ! cycle through the phase properties to find the name and species
256 ! and load the remaining data into the phase property set
257 next => null()
258 call json%get_child(j_obj, child)
259 num_spec = 0
260 do while (associated(child))
261 call json%info(child, name=key, var_type=var_type)
262 ! phase name
263 if (key.eq."name") then
264 if (var_type.ne.json_string) call die_msg(429142134, &
265 "Received non-string aerosol phase name.")
266 call json%get(child, unicode_str_val)
267 this%phase_name = unicode_str_val
268 ! chemical species in the phase
269 else if (key.eq."species") then
270 if (var_type.ne.json_array) call die_msg(293312378, &
271 "Received non-array list of aerosol phase species: "//&
272 to_string(var_type))
273 call json%get_child(child, species)
274 do while (associated(species))
275 num_spec = num_spec + 1
276 call json%info(species, var_type=var_type)
277 if (var_type.eq.json_object) then
278 ! handle species object
279 call json%get(species, "name", species_name)
280 call this%add(species_name)
281 else if (var_type.eq.json_string) then
282 ! string species name
283 call json%get(species, unicode_str_val)
284 str_val = unicode_str_val
285 call this%add(str_val)
286 else
287 call die_msg(391082805, "Invalid species format: must be object or string.")
288 end if
289 call json%get_next(species, next)
290 species => next
291 end do
292
293 ! load remaining properties into the phase property set
294 else if (key.ne."type") then
295 call property_set%load(json, child, .false., this%phase_name)
296 end if
297
298 call json%get_next(child, next)
299 child => next
300 end do
301
302 ! save the property set
303 if (associated(this%property_set)) then
304 call this%property_set%update(property_set, this%phase_name)
305 deallocate (property_set)
306 else
307 this%property_set => property_set
308 end if
309
310 ! cycle through the species associated with each phase and add
311 ! associated properties
312 i_spec = 0
313 ! allocate space for species property sets associated with a phase
314 allocate(this%spec_property_set(num_spec))
315 next => null()
316 call json%get_child(j_obj, child)
317 do while (associated(child))
318 call json%info(child, name=key, var_type=var_type)
319 ! chemical species in the phase
320 if (key.eq."species") then
321 call json%get_child(child, species)
322 do while (associated(species))
323 i_spec = i_spec + 1
324 allocate(spec_property_set)
325 call json%info(species, var_type=var_type)
326 if (var_type.eq.json_object) then
327 call json%get_child(species, species_child)
328 do while (associated(species_child))
329 call json%info(species_child, name=key, var_type=var_type)
330 ! load remaining properties into the species property set
331 if (key.ne."name".and.key.ne."type") then
332 call spec_property_set%load(json, species_child, .false., this%spec_name(i_spec)%string)
333 end if
334 call json%get_next(species_child, next)
335 species_child => next
336 end do
337 this%spec_property_set(i_spec)%val_ => spec_property_set
338 spec_property_set => null()
339 else if (var_type.eq.json_string) then
340 ! species given as just a string name → still give an empty set
341 this%spec_property_set(i_spec)%val_ => spec_property_set
342 spec_property_set => null()
343 else
344 deallocate(spec_property_set)
345 call die_msg(195829403, "Invalid type for Phase species")
346 end if
347 call json%get_next(species, next)
348 species => next
349 end do
350 end if
351 call json%get_next(child, next)
352 child => next
353 end do
354
355#else
356 subroutine load(this)
357
358 !> Aerosol phase data
359 class(aero_phase_data_t), intent(in) :: this
360
361 call warn_msg(236665532, "No support for input files.")
362#endif
363 end subroutine load
364
365!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
366
367 !> Initialize the aerosol phase data, validating species names.
368 subroutine initialize(this, chem_spec_data)
369
370 !> Aerosol phase data
371 class(aero_phase_data_t), intent(inout) :: this
372 !> Chemical species data
373 type(chem_spec_data_t), target, intent(in) :: chem_spec_data
374
375 type(property_t), pointer :: spec_props
376 integer(kind=i_kind) :: i_spec, i_spec_phase_type
377 character(len=:), allocatable :: key_name
378
379 ! Allocate space for the condensed data arrays
380 allocate(this%condensed_data_int(num_int_prop_+this%num_spec))
381 allocate(this%condensed_data_real(num_real_prop_+2*this%num_spec))
382
383 ! Set the number of species
384 num_state_var_ = this%num_spec
385
386 ! Find the aerosol-phase species, and save their needed properties
387 do i_spec = 1, num_state_var_
388
389 ! Get the species properties
390 call assert_msg(140971956, chem_spec_data%get_property_set( &
391 this%spec_name(i_spec)%string, spec_props), &
392 "Missing property set for species '"// &
393 this%spec_name(i_spec)%string// &
394 "' in aerosol phase '"//this%phase_name//"'")
395
396 ! Get the species type
397 call assert_msg(129442398, chem_spec_data%get_type( &
398 this%spec_name(i_spec)%string, spec_type_(i_spec)), &
399 "Missing type for species '"// &
400 this%spec_name(i_spec)%string// &
401 "' in aerosol phase '"//this%phase_name//"'")
402
403 ! Make sure the species is an aerosol-phase species
404 call assert_msg(136357145, chem_spec_data%get_phase( &
405 this%spec_name(i_spec)%string, i_spec_phase_type ), &
406 "Error getting phase for species "// &
407 this%spec_name(i_spec)%string)
408 call assert_msg(861388228, i_spec_phase_type.eq.chem_spec_aero_phase, &
409 "Trying to add non-aerosol phase species to aerosol phase "// &
410 this%phase_name//"; species: "//this%spec_name(i_spec)%string)
411
412 ! get the molecular weight and density of species
413 ! present in the phase
414 if (spec_type_(i_spec).eq.chem_spec_variable .or. &
415 spec_type_(i_spec).eq.chem_spec_constant .or. &
416 spec_type_(i_spec).eq.chem_spec_pssa) then
417
418 ! Get the molecular weight
419 key_name = "molecular weight [kg mol-1]"
420 call assert_msg(512254139, &
421 spec_props%get_real(key_name, mw_(i_spec)), &
422 "Missing molecular weight for species '"// &
423 this%spec_name(i_spec)%string// &
424 "' in aerosol phase '"//this%phase_name//"'")
425
426 ! Get the density
427 key_name = "density [kg m-3]"
428 call assert_msg(224966878, &
429 spec_props%get_real(key_name, density_(i_spec)), &
430 "Missing density for species '"// &
431 this%spec_name(i_spec)%string// &
432 "' in aerosol phase '"//this%phase_name//"'")
433
434 ! activity coefficients do not need molecular weight or density
435 else
436
437 mw_(i_spec) = 0.0
438 density_(i_spec) = 0.0
439
440 end if
441
442 end do
443
444 ! Save a pointer to the chemical species data
445 this%chem_spec_data => chem_spec_data
446
447 end subroutine initialize
448
449!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
450
451 !> Get the aerosol phase name
452 function get_name(this) result (phase_name)
453
454 !> The name of the aerosol phase
455 character(len=:), allocatable :: phase_name
456 !> Aerosol phase data
457 class(aero_phase_data_t), intent(in) :: this
458
459 phase_name = this%phase_name
460
461 end function get_name
462
463!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
464
465 !> Get the number of species in the phase
466 integer(kind=i_kind) function get_size(this) result(num_spec)
467
468 !> Aerosol phase data
469 class(aero_phase_data_t), intent(in) :: this
470
471 num_spec = this%num_spec
472
473 end function get_size
474
475!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
476
477 !> Get the number of Jacobian row elements needed during solving
478 integer(kind=i_kind) function num_jac_elem(this) result(num_elem)
479
480 !> Aerosol phase data
481 class(aero_phase_data_t), intent(in) :: this
482
483 integer :: i_spec
484
485 num_elem = 0
486 do i_spec = 1, this%num_spec
487 if (spec_type_(i_spec) == chem_spec_variable .or. &
488 spec_type_(i_spec) == chem_spec_constant .or. &
489 spec_type_(i_spec) == chem_spec_pssa) then
490 num_elem = num_elem + 1
491 end if
492 end do
493
494 end function num_jac_elem
495
496!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
497
498 !> Get the aerosol phase property set
499 function get_property_set(this) result (property_set)
500
501 !> A pointer to the aerosol phase property set
502 class(property_t), pointer :: property_set
503 !> Aerosol phase data
504 class(aero_phase_data_t), intent(in) :: this
505
506 property_set => this%property_set
507
508 end function get_property_set
509
510!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
511
512 !> Get the aerosol phase species property set
513 function get_spec_property_set(this, spec_name) result (spec_property_set)
514
515 !> A pointer to the aerosol phase property set
516 class(property_t), pointer :: spec_property_set
517 !> Species name to find properties of
518 character(len=*), intent(in) :: spec_name
519 !> Aerosol phase data
520 class(aero_phase_data_t), intent(in) :: this
521
522 integer(i_kind) :: i_spec
523
524 i_spec = this%find(spec_name)
525 spec_property_set => this%spec_property_set(i_spec)%val_
526
527 end function get_spec_property_set
528
529!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
530
531 !> Get an aerosol phase species name
532 function get_species_names(this) result (spec_names)
533
534 !> Names of species in this phase
535 type(string_t), allocatable :: spec_names(:)
536 !> Aerosol phase data
537 class(aero_phase_data_t), intent(in) :: this
538
539 integer(kind=i_kind) :: i_spec
540
541 allocate(spec_names(this%num_spec))
542 do i_spec = 1, this%num_spec
543 spec_names(i_spec)%string = this%spec_name(i_spec)%string
544 end do
545
546 end function get_species_names
547
548!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
549
550 !> Get an aerosol phase species type
551 function get_species_type(this, spec_name) result (spec_type)
552
553 !> The type of a species in this phase
554 integer(kind=i_kind) :: spec_type
555 !> Aerosol phase data
556 class(aero_phase_data_t), intent(in) :: this
557 !> Name of the species
558 character(len=*), intent(in) :: spec_name
559
560 call assert_msg(163269315, this%find(spec_name).gt.0, &
561 "Species '"//spec_name//"' is not in aerosol phase '"// &
562 this%phase_name//"'.")
563 call assert(255786656, this%chem_spec_data%get_type(spec_name, spec_type))
564
565 end function get_species_type
566
567!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
568
569 !> Determine the size of a binary required to pack the aerosol
570 !! representation data
571 integer(kind=i_kind) function pack_size(this, comm)
572
573 !> Aerosol representation data
574 class(aero_phase_data_t), intent(in) :: this
575 !> MPI communicator
576 integer, intent(in) :: comm
577
578 pack_size = &
579 camp_mpi_pack_size_real_array(this%condensed_data_real, comm) + &
580 camp_mpi_pack_size_integer_array(this%condensed_data_int, comm)
581
582 end function pack_size
583
584!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
585
586 !> Pack the given value to the buffer, advancing position
587 subroutine bin_pack(this, buffer, pos, comm)
588
589 !> Aerosol representation data
590 class(aero_phase_data_t), intent(in) :: this
591 !> Memory buffer
592 character, intent(inout) :: buffer(:)
593 !> Current buffer position
594 integer, intent(inout) :: pos
595 !> MPI communicator
596 integer, intent(in) :: comm
597
598#ifdef CAMP_USE_MPI
599 integer :: prev_position
600
601 prev_position = pos
602 call camp_mpi_pack_real_array(buffer, pos, this%condensed_data_real, comm)
603 call camp_mpi_pack_integer_array(buffer, pos, this%condensed_data_int,comm)
604 call assert(561436372, &
605 pos - prev_position <= this%pack_size(comm))
606#endif
607
608 end subroutine bin_pack
609
610!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
611
612 !> Unpack the given value from the buffer, advancing position
613 subroutine bin_unpack(this, buffer, pos, comm)
614
615 !> Aerosol representation data
616 class(aero_phase_data_t), intent(out) :: this
617 !> Memory buffer
618 character, intent(inout) :: buffer(:)
619 !> Current buffer position
620 integer, intent(inout) :: pos
621 !> MPI communicator
622 integer, intent(in) :: comm
623
624#ifdef CAMP_USE_MPI
625 integer :: prev_position
626
627 prev_position = pos
628 call camp_mpi_unpack_real_array(buffer, pos, this%condensed_data_real,comm)
629 call camp_mpi_unpack_integer_array(buffer, pos, this%condensed_data_int, &
630 comm)
631 call assert(219217030, &
632 pos - prev_position <= this%pack_size(comm))
633#endif
634
635 end subroutine bin_unpack
636
637!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
638
639 !> Print out the aerosol phase data
640 subroutine do_print(this, file_unit)
641
642 !> Aerosol phase data
643 class(aero_phase_data_t), intent(in) :: this
644 !> File unit for output
645 integer(kind=i_kind), optional :: file_unit
646
647 integer(kind=i_kind) :: f_unit
648 integer(kind=i_kind) :: i_spec
649
650 f_unit = 6
651
652 if (present(file_unit)) f_unit = file_unit
653 write(f_unit,*) "Aerosol phase: ", this%phase_name
654 write(f_unit,*) "Number of species: ", this%num_spec
655 write(f_unit,*) "Species: ["
656 do i_spec = 1, this%num_spec
657 write(*,*) this%spec_name(i_spec)%string
658 end do
659 write(f_unit,*) "]"
660 call this%property_set%print(f_unit)
661 write(f_unit,*) "End aerosol phase: ", this%phase_name
662
663 end subroutine do_print
664
665!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
666
667 !> Finalize the aerosol phase data
668 subroutine finalize(this)
669
670 !> Aerosol phase data
671 type(aero_phase_data_t), intent(inout) :: this
672 integer(kind=i_kind) :: i
673
674 if (allocated(this%phase_name)) deallocate(this%phase_name)
675 if (associated(this%spec_name)) deallocate(this%spec_name)
676 if (associated(this%property_set)) deallocate(this%property_set)
677 if (allocated(this%condensed_data_real)) &
678 deallocate(this%condensed_data_real)
679 if (allocated(this%condensed_data_int)) &
680 deallocate(this%condensed_data_int)
681
682 if (allocated(this%spec_property_set)) then
683 do i = 1, size(this%spec_property_set)
684 if (associated(this%spec_property_set(i)%val_)) then
685 deallocate(this%spec_property_set(i)%val_)
686 end if
687 end do
688 deallocate(this%spec_property_set)
689 end if
690
691 end subroutine finalize
692
693!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
694
695 !> Finalize the aerosol phase data
696 subroutine finalize_array(this)
697
698 !> Aerosol phase data
699 type(aero_phase_data_t), intent(inout) :: this(:)
700
701 integer(kind=i_kind) :: i_phase
702
703 do i_phase = 1, size(this)
704 call finalize(this(i_phase))
705 end do
706
707 end subroutine finalize_array
708
709!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
710
711 !> Ensure there is enough room in the species dataset to add a specified
712 !! number of species
713 subroutine ensure_size(this, num_spec)
714
715 !> Aerosol phase data
716 class(aero_phase_data_t), intent(inout) :: this
717 !> Number of new species to ensure space for
718 integer(kind=i_kind), intent(in) :: num_spec
719
720 integer :: new_size
721 type(string_t), pointer :: new_name(:)
722
723 if (size(this%spec_name) .ge. this%num_spec + num_spec) return
724 new_size = this%num_spec + num_spec + realloc_inc
725 allocate(new_name(new_size))
726 new_name(1:this%num_spec) = this%spec_name(1:this%num_spec)
727 deallocate(this%spec_name)
728 this%spec_name => new_name
729
730 end subroutine ensure_size
731
732!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
733
734 !> Add a new chemical species to the phase
735 subroutine add(this, spec_name)
736
737 !> Aerosol phase data
738 class(aero_phase_data_t), intent(inout) :: this
739 !> Name of the species to add
740 character(len=*), intent(in) :: spec_name
741
742 integer(kind=i_kind) :: i_spec
743
744 i_spec = this%find(spec_name)
745 if (i_spec.ne.0) then
746 call warn_msg(980242449, "Species "//spec_name//&
747 " added more than once to phase "//this%name())
748 return
749 end if
750 call this%ensure_size(1)
751 this%num_spec = this%num_spec + 1
752 this%spec_name(this%num_spec)%string = spec_name
753
754 end subroutine add
755
756!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
757
758 !> Get the index of an aerosol-phase species by name. Return 0 if the
759 !! species is not found
760 integer(kind=i_kind) function find(this, spec_name) &
761 result(spec_id)
762
763 !> Aerosol phase data
764 class(aero_phase_data_t), intent(in) :: this
765 !> Species name
766 character(len=*), intent(in) :: spec_name
767
768 integer(kind=i_kind) :: i_spec
769
770 spec_id = 0
771 do i_spec = 1, this%num_spec
772 if (this%spec_name(i_spec)%string .eq. spec_name) then
773 spec_id = i_spec
774 return
775 end if
776 end do
777
778 end function find
779
780!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
781
782 !> Dereference a pointer to aerosol phase data
783 elemental subroutine dereference(this)
784
785 !> Pointer to aerosol phase data
786 class(aero_phase_data_ptr), intent(inout) :: this
787
788 this%val => null()
789
790 end subroutine dereference
791
792!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
793
794 !> Finalize a pointer to aerosol phase data
795 subroutine ptr_finalize(this)
796
797 !> Pointer to aerosol phase data
798 type(aero_phase_data_ptr), intent(inout) :: this
799
800 if (associated(this%val)) deallocate(this%val)
801
802 end subroutine ptr_finalize
803
804!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
805
806 !> Finalize an array of pointers to aerosol phase data
807 subroutine ptr_finalize_array(this)
808
809 !> Array of pointers to aerosol phase data
810 type(aero_phase_data_ptr), intent(inout) :: this(:)
811
812 integer(kind=i_kind) :: i_ptr
813
814 do i_ptr = 1, size(this)
815 call ptr_finalize(this(i_ptr))
816 end do
817
818 end subroutine ptr_finalize_array
819
820!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
821
822#undef NUM_STATE_VAR_
823#undef NUM_INT_PROP_
824#undef NUM_REAL_PROP_
825#undef SPEC_TYPE_
826#undef MW_
827#undef DENSITY_
828
829end 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.
subroutine finalize_array(this)
Finalize the aerosol phase data.
subroutine bin_unpack(this, buffer, pos, comm)
Unpack the given value from the buffer, advancing position.
class(property_t) function, pointer get_spec_property_set(this, spec_name)
Get the aerosol phase species property set.
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 ptr_finalize_array(this)
Finalize an array of pointers to aerosol phase data.
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.
subroutine finalize(this)
Finalize the aerosol phase data.
integer(kind=i_kind) function pack_size(this, comm)
Determine the size of a binary required to pack the aerosol representation data.
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:53