CAMP 1.0.0
Chemistry Across Multiple Phases
aero_rep_modal_binned_mass.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_rep_modal_binned_mass module.
7
8!> \page camp_aero_rep_modal_binned_mass CAMP: Modal/Binned Mass Aerosol Representation
9!!
10!! The modal/binned mass aerosol representation includes a set of sections/bins
11!! that are made up of one or more \ref camp_aero_phase "aerosol phases." The
12!! \c json object for this \ref camp_aero_rep "aerosol representation" has the
13!! following format:
14!! \code{.json}
15!! { "camp-data" : [
16!! {
17!! "name" : "my modal/binned aero rep",
18!! "type" : "AERO_REP_MODAL_BINNED_MASS",
19!! "modes/bins" :
20!! {
21!! "dust" :
22!! {
23!! "type" : "BINNED",
24!! "phases" : [ "insoluble", "organic", "aqueous" ],
25!! "bins" : 8,
26!! "minimum diameter [m]" : 0.8e-9,
27!! "maximum diameter [m]" : 1.0e-6,
28!! "scale" : "LOG"
29!! },
30!! "depeche" :
31!! {
32!! "type" : "MODAL",
33!! "phases" : [ "moody", "listless" ],
34!! "shape" : "LOG_NORMAL",
35!! }
36!! }
37!! },
38!! ...
39!! ]}
40!! \endcode
41!! The key-value pair \b type is required and must be
42!! \b AERO_REP_MODAL_BINNED_MASS. The key-value pair \b modes/bins is also
43!! required and must contain a set of at least one uniquely named mode or
44!! bin-set key-value pair whose value(s) specify a \b type that must be either
45!! \b MODAL or \b BINNED and an array of \b phases that correspond to existing
46!! \ref camp_aero_phase "aerosol phase" objects. Each phase will be present
47!! once within a mode or once within each bin in a bin-set.
48!!
49!! Modes must also specify a distribution \b shape which must be \b LOG_NORMAL
50!! (the available shapes may be expanded in the future). Log-normal sections
51!! have \b geometric \b mean \b diameter (m) and a \b geometric
52!! \b standard \b deviation (unitless) properties that must be set using an
53!! \c aero_rep_update_data_modal_binned_mass_GMD_t or and
54!! \c aero_rep_update_data_modal_binned_mass_GSD_t object and will be used along
55!! with the mass
56!! concentration of species in each phase and their densities to calculate a
57!! lognormal distribution for each mode at runtime.
58!!
59!! Bin sets must specify the number of \b bins, a \b minimum \b diameter (m),
60!! a \b maximum \b diameter (m) and a \b scale, which must be \b LOG or
61!! \b LINEAR. The number concentration will be calculated at run-time based on
62!! the total mass of each bin, the species densities and the diameter of
63!! particles in that bin.
64!!
65!! The GMD and GSD for each mode must be set from an external model using
66!! \c camp_aero_rep_modal_binned_mass::aero_rep_update_data_modal_binned_mass_GMD_t
67!! and
68!! \c camp_aero_rep_modal_binned_mass::aero_rep_update_data_modal_binned_mass_GSD_t
69!! objects.
70
71!> The abstract aero_rep_modal_binned_mass_t structure and associated subroutines.
73
78 use camp_mpi
80 use camp_util, only: dp, i_kind, &
83
84 use iso_c_binding
85
86 implicit none
87 private
88
89#define BINNED 1
90#define MODAL 2
91
92#define NUM_SECTION_ this%condensed_data_int(1)
93#define INT_DATA_SIZE_ this%condensed_data_int(2)
94#define REAL_DATA_SIZE_ this%condensed_data_int(3)
95#define AERO_REP_ID_ this%condensed_data_int(4)
96#define NUM_INT_PROP_ 4
97#define NUM_REAL_PROP_ 0
98#define NUM_ENV_PARAM_ 0
99#define MODE_INT_PROP_LOC_(x) this%condensed_data_int(NUM_INT_PROP_+x)
100#define MODE_REAL_PROP_LOC_(x) this%condensed_data_int(NUM_INT_PROP_+NUM_SECTION_+x)
101#define SECTION_TYPE_(x) this%condensed_data_int(MODE_INT_PROP_LOC_(x))
102
103! For modes, NUM_BINS_ = 1
104#define NUM_BINS_(x) this%condensed_data_int(MODE_INT_PROP_LOC_(x)+1)
105
106! Number of aerosol phases in this mode/bin set
107#define NUM_PHASE_(x) this%condensed_data_int(MODE_INT_PROP_LOC_(x)+2)
108
109! Phase state and model data ids
110#define PHASE_STATE_ID_(x,y,b) this%condensed_data_int(MODE_INT_PROP_LOC_(x)+2+(b-1)*NUM_PHASE_(x)+y)
111#define PHASE_MODEL_DATA_ID_(x,y,b) this%condensed_data_int(MODE_INT_PROP_LOC_(x)+2+NUM_BINS_(x)*NUM_PHASE_(x)+(b-1)*NUM_PHASE_(x)+y)
112
113! Number of Jacobian elements in a phase
114#define PHASE_NUM_JAC_ELEM_(x,y,b) this%condensed_data_int(MODE_INT_PROP_LOC_(x)+2+2*NUM_BINS_(x)*NUM_PHASE_(x)+(b-1)*NUM_PHASE_(x)+y)
115
116! Bin diameter (bins only)
117#define BIN_DP_(x,b) this%condensed_data_real(MODE_REAL_PROP_LOC_(x)+(b-1)*3)
118
119! Real-time number concetration - used for modes and bins - for modes, b=1
120#define NUMBER_CONC_(x,b) this%condensed_data_real(MODE_REAL_PROP_LOC_(x)+(b-1)*3+1)
121
122! Real-time effective radius - for modes, b=1
123#define EFFECTIVE_RADIUS_(x,b) this%condensed_data_real(MODE_REAL_PROP_LOC_(x)+(b-1)*3+2)
124
125! Real-time aerosol phase mass - used for modes and bins - for modes, b=1
126#define PHASE_MASS_(x,y,b) this%condensed_data_real(MODE_REAL_PROP_LOC_(x)+3*NUM_BINS_(x)+(b-1)*NUM_PHASE_(x)+y-1)
127
128! Real-time aerosol phase average MW - used for modes and bins - for modes, b=0
129#define PHASE_AVG_MW_(x,y,b) this%condensed_data_real(MODE_REAL_PROP_LOC_(x)+(3+NUM_PHASE_(x))*NUM_BINS_(x)+(b-1)*NUM_PHASE_(x)+y-1)
130
131 ! Update types (These must match values in aero_rep_modal_binned_mass.c)
132 integer(kind=i_kind), parameter, public :: update_gmd = 0
133 integer(kind=i_kind), parameter, public :: update_gsd = 1
134
138
139 !> Modal mass aerosol representation
140 !!
141 !! Time-invariant data related to a modal/binned mass aerosol representation.
143 !> Mode names (only used during initialization)
144 type(string_t), allocatable :: section_name(:)
145 !> Phase state id (only used during initialization)
146 integer(kind=i_kind), allocatable :: phase_state_id(:)
147 contains
148 !> Initialize the aerosol representation data, validating component data and
149 !! loading any required information from the \c
150 !! aero_rep_data_t::property_set. This routine should be called once for
151 !! each aerosol representation at the beginning of a model run after all
152 !! the input files have been read in. It ensures all data required during
153 !! the model run are included in the condensed data arrays.
154 procedure :: initialize
155 !> Initialize an update data GSD object
156 procedure :: update_data_initialize_gsd => update_data_init_gsd
157 !> Initialize an update data GMD object
158 procedure :: update_data_initialize_gmd => update_data_init_gmd
159 !> Get an id for a mode or bin in the aerosol representation by name for
160 !! use with updates from external modules
161 procedure :: get_section_id
162 !> Get the size of the section of the
163 !! \c camp_camp_state::camp_state_t::state_var array required for this
164 !! aerosol representation.
165 !!
166 !! For a modal/binned mass representation, the size will correspond to the
167 !! the sum of the sizes of a single instance of each aerosol phase
168 !! provided to \c aero_rep_modal_binned_mass::initialize()
169 procedure :: size => get_size
170 !> Get a list of unique names for each element on the
171 !! \c camp_camp_state::camp_state_t::state_var array for this aerosol
172 !! representation. The list may be restricted to a particular phase and/or
173 !! aerosol species by including the phase_name and spec_name arguments.
174 !!
175 !! For a modal/binned mass representation, the unique names for bins are:
176 !! - "bin name.bin #.phase name.species name"
177 !!
178 !! ... and for modes are:
179 !! - "mode name.phase name.species name"
180 procedure :: unique_names
181 !> Get a species id on the \c camp_camp_state::camp_state_t::state_var
182 !! array by its unique name. These are unique ids for each element on the
183 !! state array for this \ref camp_aero_rep "aerosol representation" and
184 !! are numbered:
185 !!
186 !! \f[x_u \in x_f ... (x_f+n-1)\f]
187 !!
188 !! where \f$x_u\f$ is the id of the element corresponding to the species
189 !! with unique name \f$u\f$ on the \c
190 !! camp_camp_state::camp_state_t::state_var array, \f$x_f\f$ is the index
191 !! of the first element for this aerosol representation on the state array
192 !! and \f$n\f$ is the total number of variables on the state array from
193 !! this aerosol representation.
194 procedure :: spec_state_id
195 !> Get the non-unique name of a species by its unique name
196 procedure :: spec_name
197 !> Get the number of instances of an aerosol phase in this representation
199 !> Get the number of Jacobian elements used in calculations of aerosol mass,
200 !! volume, number, etc. for a particular phase
201 procedure :: num_jac_elem
202 !> Finalize the aerosol representation
203 final :: finalize
204
206
207 ! Constructor for aero_rep_modal_binned_mass_t
209 procedure :: constructor
210 end interface aero_rep_modal_binned_mass_t
211
212 !> Update GMD object
213 type, extends(aero_rep_update_data_t) :: &
215 private
216 !> Flag indicating whether the update data has been allocated
217 logical :: is_malloced = .false.
218 !> Unique id for finding aerosol representations during initialization
219 integer(kind=i_kind) :: aero_rep_unique_id = 0
220 contains
221 !> Update the GMD
222 procedure :: set_gmd => update_data_set_gmd
223 !> Determine the pack size of the local update data
225 !> Pack the local update data to a binary
227 !> Unpack the local update data from a binary
229 !> Finalize the GMD update data
232
233 !> Update GSD object
234 type, extends(aero_rep_update_data_t) :: &
236 private
237 !> Flag indicating whether the update data has been allocated
238 logical :: is_malloced = .false.
239 !> Unique id for finding aerosol representations during initialization
240 integer(kind=i_kind) :: aero_rep_unique_id = 0
241 contains
242 !> Update the GSD
243 procedure :: set_gsd => update_data_set_gsd
244 !> Determine the pack size of the local update data
246 !> Pack the local update data to a binary
248 !> Unpack the local update data from a binary
250 !> Finalize the GSD update data
253
254 !> Interface to c aerosol representation functions
255 interface
256
257 !> Allocate space for a GMD update object
259 result(update_data) bind (c)
260 use iso_c_binding
261 !> Allocated update data object
262 type(c_ptr) :: update_data
264
265 !> Set a new mode GMD
267 aero_rep_unique_id, section_id, gmd) bind (c)
268 use iso_c_binding
269 !> Update data
270 type(c_ptr), value :: update_data
271 !> Aerosol representation unique id
272 integer(kind=c_int), value :: aero_rep_unique_id
273 !> Section id from
274 !! camp_aero_rep_modal_binned_mass::aero_rep_modal_binned_mass_t::get_section_id
275 integer(kind=c_int), value :: section_id
276 !> New GMD (m)
277 real(kind=c_double), value :: gmd
279
280 !> Allocate space for a GSD update object
282 result(update_data) bind (c)
283 use iso_c_binding
284 !> Allocated update data object
285 type(c_ptr) :: update_data
287
288 !> Set a new mode GSD
290 aero_rep_unique_id, section_id, gsd) bind (c)
291 use iso_c_binding
292 !> Update data
293 type(c_ptr), value :: update_data
294 !> Aerosol representation unique id
295 integer(kind=c_int), value :: aero_rep_unique_id
296 !> Section id from
297 !! camp_aero_rep_modal_binned_mass::aero_rep_modal_binned_mass_t::get_section_id
298 integer(kind=c_int), value :: section_id
299 !> New GSD (m)
300 real(kind=c_double), value :: gsd
302
303 !> Free an update data object
304 pure subroutine aero_rep_free_update_data(update_data) bind (c)
305 use iso_c_binding
306 !> Update data
307 type(c_ptr), value, intent(in) :: update_data
308 end subroutine aero_rep_free_update_data
309
310 end interface
311
312contains
313
314!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
315
316 !> Constructor for aero_rep_modal_binned_mass_t
317 function constructor() result (new_obj)
318
319 !> New aerosol representation
320 type(aero_rep_modal_binned_mass_t), pointer :: new_obj
321
322 allocate(new_obj)
323
324 end function constructor
325
326!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
327
328 !> Initialize the aerosol representation data, validating component data and
329 !! loading any required information from the \c
330 !! aero_rep_data_t::property_set. This routine should be called once for
331 !! each aerosol representation at the beginning of a model run after all
332 !! the input files have been read in. It ensures all data required during
333 !! the model run are included in the condensed data arrays.
334 subroutine initialize(this, aero_phase_set, spec_state_id)
335
336 !> Aerosol representation data
337 class(aero_rep_modal_binned_mass_t), intent(inout) :: this
338 !> The set of aerosol phases
339 type(aero_phase_data_ptr), pointer, intent(in) :: aero_phase_set(:)
340 !> Beginning state id for this aerosol representationin the model species
341 !! state array
342 integer(kind=i_kind), intent(in) :: spec_state_id
343
344 type(property_t), pointer :: sections, section, phases
345 integer(kind=i_kind) :: i_section, i_phase, j_phase, k_phase, &
346 i_bin
347 integer(kind=i_kind) :: curr_spec_state_id
348 integer(kind=i_kind) :: num_phase, num_bin
349 integer(kind=i_kind) :: n_int_param, n_float_param
350 character(len=:), allocatable :: key_name, phase_name, sect_type, str_val
351 real(kind=dp) :: min_dp, max_dp, d_log_dp
352
353 ! Determine the size of the condensed data arrays
354 n_int_param = num_int_prop_
355 n_float_param = num_real_prop_
356
357 ! Get the set of sections/bin-sets
358 key_name = "modes/bins"
359 call assert_msg(877855909, &
360 this%property_set%get_property_t(key_name, sections), &
361 "Missing sections/bins for modal/binned mass aerosol "// &
362 "representation '"//this%rep_name//"'")
363 call assert_msg(894962494, sections%size().gt.0, "No sections or bins "// &
364 "specified for modal/binned mass aerosol representation '"// &
365 this%rep_name//"'")
366
367 ! Allocate space for the mode/bin names
368 allocate(this%section_name(sections%size()))
369
370 ! Loop through the sections, adding names and counting the spaces needed
371 ! on the condensed data arrays, and counting the total phases instances
372 num_phase = 0
373 call sections%iter_reset()
374 do i_section = 1, sections%size()
375
376 ! Get the mode/bin name
377 call assert(867378489, sections%get_key(key_name))
378 call assert_msg(234513113, len(key_name).gt.0, "Missing mode/bin "// &
379 "name in modal/binned mass aerosol representation '"// &
380 this%rep_name//"'")
381 this%section_name(i_section)%string = key_name
382
383 ! Get the mode/bin properties
384 call assert_msg(517138327, sections%get_property_t(val=section), &
385 "Invalid structure for mode/bin '"// &
386 this%section_name(i_section)%string// &
387 "' in modal/binned mass aerosol representation '"// &
388 this%rep_name//"'")
389
390 ! Get the section type
391 key_name = "type"
392 call assert_msg(742404898, section%get_string(key_name, sect_type), &
393 "Missing mode/bin type in mode/bin '"// &
394 this%section_name(i_section)%string// &
395 "' in modal/binned mass aerosol representation '"// &
396 this%rep_name//"'")
397 call assert_msg(618556995, &
398 sect_type.eq."MODAL".or.sect_type.eq."BINNED", &
399 "Invalid mode/bin type '"//sect_type//"' in mode/bin '"// &
400 this%section_name(i_section)%string// &
401 "' in modal/binned mass aerosol representation '"// &
402 this%rep_name//"'")
403
404 ! Get the number of bins (or set to 1 for a mode)
405 num_bin = 1
406 if (sect_type.eq."BINNED") then
407
408 key_name = "bins"
409 call assert_msg(824494286, section%get_int(key_name, num_bin), &
410 "Missing number of bins in bin '"// &
411 this%section_name(i_section)%string// &
412 "' in modal/binned mass aerosol representation '"// &
413 this%rep_name//"'")
414 end if
415
416 ! Add space for the mode/bin type, number of bins, and phase count
417 ! and parameter locations
418 n_int_param = n_int_param + 5
419
420 ! Add space for the bin diameter, number concentration, and
421 ! effective radius
422 n_float_param = n_float_param + 3*num_bin
423
424 ! Get the set of phases
425 key_name = "phases"
426 call assert_msg(815518058, section%get_property_t(key_name, phases), &
427 "Missing phases for mode '"// &
428 this%section_name(i_section)%string// &
429 "' in modal/binned mass aerosol representation '"// &
430 this%rep_name//"'")
431
432 ! Add the phases to the counter
433 call assert_msg(772593427, phases%size().gt.0, &
434 "No phases specified for mode '"// &
435 this%section_name(i_section)%string// &
436 "' in modal/binned mass aerosol representation '"// &
437 this%rep_name//"'")
438 num_phase = num_phase + phases%size() * num_bin
439
440 ! Loop through the phases and make sure they exist
441 call phases%iter_reset()
442 do i_phase = 1, phases%size()
443
444 ! Get the phase name
445 call assert_msg(393427582, phases%get_string(val=phase_name), &
446 "Non-string phase name for mode '"// &
447 this%section_name(i_section)%string// &
448 "' in modal/binned mass aerosol representation '"// &
449 this%rep_name//"'")
450
451 ! Find the aerosol phase and add space for its variables
452 do j_phase = 1, size(aero_phase_set)
453 if (phase_name.eq.aero_phase_set(j_phase)%val%name()) then
454
455 ! Add space for the phase state, model data ids, and number
456 ! of Jacobian elements
457 n_int_param = n_int_param + 3 * num_bin
458
459 ! Add space for total aerosol phase mass and average MW,
460 n_float_param = n_float_param + 2 * num_bin
461
462 exit
463 else if (j_phase.eq.size(aero_phase_set)) then
464 call die_msg(652391420, "Non-existant aerosol phase '"// &
465 phase_name//"' specified for mode '"// &
466 this%section_name(i_section)%string// &
467 "' in modal/binned mass aerosol representation '"// &
468 this%rep_name//"'")
469 end if
470 end do
471
472 call phases%iter_next()
473 end do
474
475 call sections%iter_next()
476 end do
477
478 ! Allocate space for the aerosol phases and species state ids
479 allocate(this%aero_phase(num_phase))
480 allocate(this%aero_phase_is_at_surface(num_phase))
481 allocate(this%phase_state_id(size(this%aero_phase)))
482
483 ! Allocate condensed data arrays
484 allocate(this%condensed_data_int(n_int_param))
485 allocate(this%condensed_data_real(n_float_param))
486 this%condensed_data_int(:) = int(0, kind=i_kind)
487 this%condensed_data_real(:) = real(0.0, kind=dp)
488 int_data_size_ = n_int_param
489 real_data_size_ = n_float_param
490
491 ! Set the number of sections
492 num_section_ = sections%size()
493
494 ! Save space for the environment-dependent parameters (GMD and GSD)
495 this%num_env_params = 2 * num_section_
496
497 ! Loop through the sections, adding names and distribution parameters and
498 ! counting the phases in each section
499 i_phase = 1
500 curr_spec_state_id = spec_state_id
501 n_int_param = num_int_prop_+2*num_section_+1
502 n_float_param = num_real_prop_+1
503 call sections%iter_reset()
504 do i_section = 1, num_section_
505
506 ! Set the data locations for this mode
507 mode_int_prop_loc_(i_section) = n_int_param
508 mode_real_prop_loc_(i_section) = n_float_param
509
510 ! Get the mode/bin properties
511 call assert(394743663, sections%get_property_t(val=section))
512
513 ! Get the mode/bin type
514 key_name = "type"
515 call assert(667058653, section%get_string(key_name, sect_type))
516 if (sect_type.eq."MODAL") then
517 section_type_(i_section) = modal
518 else if (sect_type.eq."BINNED") then
519 section_type_(i_section) = binned
520 else
521 call die_msg(256924433, "Internal error")
522 end if
523
524 ! Get the number of bins (or set to 1 for a mode)
525 num_bins_(i_section) = 1
526 if (section_type_(i_section).eq.binned) then
527 key_name = "bins"
528 call assert(315215287, section%get_int(key_name, num_bins_(i_section)))
529 end if
530
531 ! Get mode parameters
532 if (section_type_(i_section).eq.modal) then
533
534 ! Currently no mode parameters
535
536 effective_radius_(i_section,1) = -9999.9
537 number_conc_(i_section,1) = -9999.9
538
539 ! Get bin parameters
540 else if (section_type_(i_section).eq.binned) then
541
542 ! Get the minimum diameter (m)
543 key_name = "minimum diameter [m]"
544 call assert_msg(548762180, section%get_real(key_name, min_dp), &
545 "Missing minimum diameter for bin '"// &
546 this%section_name(i_section)%string// &
547 "' in modal/binned mass aerosol representation '"// &
548 this%rep_name//"'")
549
550 ! Get the maximum diameter (m)
551 key_name = "maximum diameter [m]"
552 call assert_msg(288632226, section%get_real(key_name, max_dp), &
553 "Missing maximum diameter for bin '"// &
554 this%section_name(i_section)%string// &
555 "' in modal/binned mass aerosol representation '"// &
556 this%rep_name//"'")
557
558 ! Get the scale
559 key_name = "scale"
560 call assert_msg(404761639, section%get_string(key_name, str_val), &
561 "Missing bin scale for bin '"// &
562 this%section_name(i_section)%string// &
563 "' in modal/binned mass aerosol representation '"// &
564 this%rep_name//"'")
565
566 ! Assign the bin diameters
567 if (str_val.eq."LINEAR") then
568 do i_bin = 1, num_bins_(i_section)
569 bin_dp_(i_section,i_bin) = min_dp + &
570 (i_bin-1) * (max_dp-min_dp)/(num_bins_(i_section)-1)
571 end do
572 else if (str_val.eq."LOG") then
573 d_log_dp = (log10(max_dp)-log10(min_dp))/(num_bins_(i_section)-1)
574 do i_bin = 1, num_bins_(i_section)
575 bin_dp_(i_section,i_bin) = 10.0d0**( log10(min_dp) + &
576 (i_bin-1) * d_log_dp )
577 end do
578 else
579 call die_msg(236797392, "Invalid scale specified for bin '"// &
580 this%section_name(i_section)%string// &
581 "' in modal/binned mass aerosol representation '"// &
582 this%rep_name//"'")
583 end if
584 do i_bin = 1, num_bins_(i_section)
585 ! Set the effective radius
586 effective_radius_(i_section,i_bin) = bin_dp_(i_section,i_bin) / 2.0
587 number_conc_(i_section,i_bin) = -9999.9
588 end do
589 end if
590
591 ! Get the set of phases
592 key_name = "phases"
593 call assert(712411046, section%get_property_t(key_name, phases))
594
595 ! Save the number of phases
596 num_phase_(i_section) = phases%size()
597
598 ! Add space for the mode/bin type, number of bins, and number of phases
599 n_int_param = n_int_param + 3
600
601 ! Add space for bin diameter, number concentration and effective radius
602 n_float_param = n_float_param + 3 * num_bins_(i_section)
603
604 ! Loop through the phase names, look them up, and add them to the list
605 call phases%iter_reset()
606 do j_phase = 1, phases%size()
607
608 ! Get the phase name
609 call assert(775801035, phases%get_string(val=phase_name))
610
611 ! Find the aerosol phase and add it to the list
612 do k_phase = 1, size(aero_phase_set)
613 if (phase_name.eq.aero_phase_set(k_phase)%val%name()) then
614
615 ! Loop through the bins
616 do i_bin = 1, num_bins_(i_section)
617
618 ! Add the aerosol phase to the list
619 this%aero_phase(i_phase) = aero_phase_set(k_phase)
620
621 ! No species exist at surface
622 this%aero_phase_is_at_surface(i_phase) = .true.
623
624 ! Save the starting id for this phase on the state array
625 this%phase_state_id(i_phase) = curr_spec_state_id
626 phase_state_id_(i_section, j_phase, i_bin) = curr_spec_state_id
627
628 ! Increment the state id by the size of the phase
629 curr_spec_state_id = curr_spec_state_id + &
630 aero_phase_set(k_phase)%val%size()
631
632 ! Save the phase model data id
633 phase_model_data_id_(i_section, j_phase, i_bin) = k_phase
634
635 i_phase = i_phase + 1
636 end do
637
638 ! Add space for aerosol phase state, model data ids, and
639 ! number of Jacobian elements
640 n_int_param = n_int_param + 3*num_bins_(i_section)
641
642 ! Add space for aerosol phase mass and average MW
643 n_float_param = n_float_param + 2*num_bins_(i_section)
644
645 exit
646 else if (k_phase.eq.size(aero_phase_set)) then
647 call die_msg(652391420, "Internal error.")
648 end if
649 end do
650
651 call phases%iter_next()
652 end do
653
654 call sections%iter_next()
655 end do
656
657 ! Initialize the aerosol representation id
658 aero_rep_id_ = -1
659
660 ! Check the data sizes
661 call assert(951534966, i_phase-1.eq.num_phase)
662 call assert(951534966, n_int_param.eq.int_data_size_+1)
663 call assert(325387136, n_float_param.eq.real_data_size_+1)
664
665 end subroutine initialize
666
667!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
668
669 !> Get an id for a mode or bin by name for use with updates from external
670 !! modules
671 function get_section_id(this, section_name, section_id) result (found)
672
673 !> Flag indicating whether the mode/bin was found
674 logical :: found
675 !> Aerosol representation
676 class(aero_rep_modal_binned_mass_t), intent(in) :: this
677 !> Section name
678 character(len=*), intent(in) :: section_name
679 !> Section id
680 integer(kind=i_kind), intent(out) :: section_id
681
682 integer(kind=i_kind) :: i_section
683
684 call assert_msg(194186171, len(trim(section_name)).gt.0, &
685 "Trying to get section id of unnamed aerosol "// &
686 "representation.")
687
688 found = .false.
689 do i_section = 1, size(this%section_name)
690 if (this%section_name(i_section)%string.eq.trim(section_name)) then
691 found = .true.
692 section_id = i_section
693 return
694 end if
695 end do
696
697 end function get_section_id
698
699!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
700
701 !> Get the size of the section of the
702 !! \c camp_camp_state::camp_state_t::state_var array required for this
703 !! aerosol representation.
704 !!
705 !! For a modal/binned mass representation, the size will correspond to the
706 !! the sum of the sizes of a single instance of each aerosol phase
707 !! provided to \c aero_rep_modal_binned_mass::initialize()
708 function get_size(this) result (state_size)
709
710 !> Size on the state array
711 integer(kind=i_kind) :: state_size
712 !> Aerosol representation data
713 class(aero_rep_modal_binned_mass_t), intent(in) :: this
714
715 integer(kind=i_kind) :: i_phase
716
717 ! Get the total number of species across all phases
718 state_size = 0
719 do i_phase = 1, size(this%aero_phase)
720 state_size = state_size + this%aero_phase(i_phase)%val%size()
721 end do
722
723 end function get_size
724
725!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
726
727 !> Get a list of unique names for each element on the
728 !! \c camp_camp_state::camp_state_t::state_var array for this aerosol
729 !! representation. The list may be restricted to a particular phase and/or
730 !! aerosol species by including the phase_name and spec_name arguments.
731 !!
732 !! For a modal/binned mass representation, the unique names for bins are:
733 !! - "bin name.bin #.phase name.species name"
734 !!
735 !! ... and for modes are:
736 !! - "mode name.phase name.species name"
737 function unique_names(this, phase_name, tracer_type, spec_name)
738
739 !> List of unique names
740 type(string_t), allocatable :: unique_names(:)
741 !> Aerosol representation data
742 class(aero_rep_modal_binned_mass_t), intent(in) :: this
743 !> Aerosol phase name
744 character(len=*), optional, intent(in) :: phase_name
745 !> Aerosol-phase species tracer type
746 integer(kind=i_kind), optional, intent(in) :: tracer_type
747 !> Aerosol-phase species name
748 character(len=*), optional, intent(in) :: spec_name
749
750 integer(kind=i_kind) :: num_spec, i_spec, j_spec, i_phase, j_phase, &
751 i_section, i_bin
752 integer(kind=i_kind) :: curr_tracer_type
753 character(len=:), allocatable :: curr_section_name, curr_phase_name, &
754 curr_bin_str
755 type(string_t), allocatable :: spec_names(:)
756
757 ! Count the number of unique names
758 num_spec = 0
759 do i_phase = 1, size(this%aero_phase)
760
761 ! Filter by phase name
762 if (present(phase_name)) then
763 curr_phase_name = this%aero_phase(i_phase)%val%name()
764 if (phase_name.ne.curr_phase_name) cycle
765 end if
766
767 ! Filter by spec name and/or tracer type
768 if (present(spec_name).or.present(tracer_type)) then
769 spec_names = this%aero_phase(i_phase)%val%get_species_names()
770 do j_spec = 1, size(spec_names)
771 curr_tracer_type = &
772 this%aero_phase(i_phase)%val%get_species_type( &
773 spec_names(j_spec)%string)
774 if (present(spec_name)) then
775 if (spec_name.ne.spec_names(j_spec)%string) cycle
776 end if
777 if (present(tracer_type)) then
778 if (tracer_type.ne.curr_tracer_type) cycle
779 end if
780 num_spec = num_spec + 1
781 end do
782 else
783 num_spec = num_spec + this%aero_phase(i_phase)%val%size()
784 end if
785
786 end do
787
788 ! Allocate space for the unique names
789 allocate(unique_names(num_spec))
790
791 ! Loop through the modes/bin sets
792 i_phase = 1
793 i_spec = 1
794 do i_section = 1, num_section_
795
796 ! Get the current section name
797 curr_section_name = this%section_name(i_section)%string
798
799 ! Loop through the phases for this mode/bin set
800 do j_phase = 1, num_phase_(i_section)
801
802 ! Set the current phase name
803 curr_phase_name = this%aero_phase(i_phase)%val%name()
804
805 ! Filter by phase name
806 if (present(phase_name)) then
807 if (phase_name.ne.curr_phase_name) then
808 i_phase = i_phase + num_bins_(i_section)
809 cycle
810 end if
811 end if
812
813 ! Get the species names in this phase
814 spec_names = this%aero_phase(i_phase)%val%get_species_names()
815
816 ! Loop through the bins (one iteration for modes)
817 do i_bin = 1, num_bins_(i_section)
818
819 ! Set the current bin label (except for single bins or mode)
820 if (num_bins_(i_section).gt.1) then
821 curr_bin_str = trim(to_string(i_bin))//"."
822 else
823 curr_bin_str = ""
824 end if
825
826 ! Add species from this phase/bin
827 num_spec = this%aero_phase(i_phase)%val%size()
828 do j_spec = 1, num_spec
829
830 ! Filter by species name
831 if (present(spec_name)) then
832 if (spec_name.ne.spec_names(j_spec)%string) cycle
833 end if
834
835 ! Filter by species tracer type
836 if (present(tracer_type)) then
837 curr_tracer_type = &
838 this%aero_phase(i_phase)%val%get_species_type( &
839 spec_names(j_spec)%string)
840 if (tracer_type.ne.curr_tracer_type) cycle
841 end if
842
843 ! Add the unique name for this species
844 unique_names(i_spec)%string = curr_section_name//"."// &
845 curr_bin_str//curr_phase_name//'.'// &
846 spec_names(j_spec)%string
847
848 i_spec = i_spec + 1
849 end do
850
851 ! Move to the next phase instance
852 i_phase = i_phase + 1
853
854 end do
855
856 deallocate(spec_names)
857
858 end do
859 end do
860
861 end function unique_names
862
863!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
864
865 !> Get a species id on the \c camp_camp_state::camp_state_t::state_var
866 !! array by its unique name. These are unique ids for each element on the
867 !! state array for this \ref camp_aero_rep "aerosol representation" and
868 !! are numbered:
869 !!
870 !! \f[x_u \in x_f ... (x_f+n-1)\f]
871 !!
872 !! where \f$x_u\f$ is the id of the element corresponding to the species
873 !! with unique name \f$u\f$ on the \c
874 !! camp_camp_state::camp_state_t::state_var array, \f$x_f\f$ is the index
875 !! of the first element for this aerosol representation on the state array
876 !! and \f$n\f$ is the total number of variables on the state array from
877 !! this aerosol representation.
878 function spec_state_id(this, unique_name) result (spec_id)
879
880 !> Species state id
881 integer(kind=i_kind) :: spec_id
882 !> Aerosol representation data
883 class(aero_rep_modal_binned_mass_t), intent(in) :: this
884 !> Unique name
885 character(len=*), intent(in) :: unique_name
886
887 type(string_t), allocatable :: unique_names(:)
888 integer(kind=i_kind) :: i_spec
889
890 spec_id = 0
891 unique_names = this%unique_names()
892 do i_spec = 1, size(unique_names)
893 if (unique_names(i_spec)%string .eq. unique_name) then
894 spec_id = this%phase_state_id(1) + i_spec - 1
895 return
896 end if
897 end do
898 call die_msg( 105414960, "Cannot find species '"//unique_name//"'" )
899
900 end function spec_state_id
901
902!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
903
904 !> Get the non-unique name of a species by its unique name
905 function spec_name(this, unique_name)
906
907 !> Chemical species name
908 character(len=:), allocatable :: spec_name
909 !> Aerosol representation data
910 class(aero_rep_modal_binned_mass_t), intent(in) :: this
911 !> Unique name of the species in this aerosol representation
912 character(len=*), intent(in) :: unique_name
913
914 ! Indices for iterators
915 integer(kind=i_kind) :: i_spec, j_spec, i_phase
916
917 ! species names in the aerosol phase
918 type(string_t), allocatable :: spec_names(:)
919
920 ! unique name list
921 type(string_t), allocatable :: unique_names(:)
922
923 unique_names = this%unique_names()
924
925 i_spec = 1
926 do i_phase = 1, size(this%aero_phase)
927 spec_names = this%aero_phase(i_phase)%val%get_species_names()
928 do j_spec = 1, this%aero_phase(i_phase)%val%size()
929 if (unique_name.eq.unique_names(i_spec)%string) then
930 spec_name = spec_names(j_spec)%string
931 end if
932 i_spec = i_spec + 1
933 end do
934 deallocate(spec_names)
935 end do
936
937 deallocate(unique_names)
938
939 end function spec_name
940
941!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
942
943 !> Get the number of instances of a specified aerosol phase.
944 function num_phase_instances(this, phase_name)
945
946 !> Number of instances of the aerosol phase
947 integer(kind=i_kind) :: num_phase_instances
948 !> Aerosol representation data
949 class(aero_rep_modal_binned_mass_t), intent(in) :: this
950 !> Aerosol phase name
951 character(len=*), intent(in) :: phase_name
952
953 integer(kind=i_kind) :: i_phase
954
956 do i_phase = 1, size(this%aero_phase)
957 if (this%aero_phase(i_phase)%val%name().eq.phase_name) then
959 end if
960 end do
961
962 end function num_phase_instances
963
964!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
965
966 !> Get the number of Jacobian elements used in calculations of aerosol mass,
967 !! volume, number, etc. for a particular phase
968 function num_jac_elem(this, phase_id)
969
970 !> Number of Jacobian elements used
971 integer(kind=i_kind) :: num_jac_elem
972 !> Aerosol respresentation data
973 class(aero_rep_modal_binned_mass_t), intent(in) :: this
974 !> Aerosol phase id
975 integer(kind=i_kind), intent(in) :: phase_id
976
977 integer(kind=i_kind) :: i_section, i_phase, i_bin, j_phase, &
978 phase_id_to_find, section_size, &
979 section_start
980
981 phase_id_to_find = phase_id
982 section_start = 1
983 do i_section = 1, num_section_
984 section_size = num_phase_(i_section) * num_bins_(i_section)
985 if( phase_id_to_find .le. section_size ) then
986 i_phase = ( phase_id_to_find - 1 ) / num_bins_(i_section) + 1
987 i_bin = mod( phase_id_to_find - 1, num_bins_(i_section) ) + 1
988 num_jac_elem = 0
989 do j_phase = section_start + i_bin - 1, &
990 section_start + i_bin - 1 + &
991 ( num_phase_(i_section) - 1 ) * num_bins_(i_section), &
992 num_bins_(i_section)
994 this%aero_phase( j_phase )%val%num_jac_elem( )
995 end do
996 return
997 end if
998 phase_id_to_find = phase_id_to_find - section_size
999 section_start = section_start + section_size
1000 end do
1001
1002 end function num_jac_elem
1003
1004!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1005
1006 !> Finalize the aerosol representation
1007 elemental subroutine finalize(this)
1008
1009 !> Aerosol representation data
1010 type(aero_rep_modal_binned_mass_t), intent(inout) :: this
1011
1012 if (allocated(this%rep_name)) deallocate(this%rep_name)
1013 if (allocated(this%aero_phase)) then
1014 ! The core will deallocate the aerosol phases
1015 call this%aero_phase(:)%dereference()
1016 deallocate(this%aero_phase)
1017 end if
1018 if (associated(this%property_set)) &
1019 deallocate(this%property_set)
1020 if (allocated(this%condensed_data_real)) &
1021 deallocate(this%condensed_data_real)
1022 if (allocated(this%condensed_data_int)) &
1023 deallocate(this%condensed_data_int)
1024
1025 end subroutine finalize
1026
1027!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1028
1029 !> Initialize a GMD update object
1030 subroutine update_data_init_gmd(this, update_data, aero_rep_type)
1031
1032 use camp_rand, only : generate_int_id
1033
1034 !> Aerosol representation to update
1035 class(aero_rep_modal_binned_mass_t), intent(inout) :: this
1036 !> Update data object
1037 class(aero_rep_update_data_modal_binned_mass_gmd_t), intent(out) :: &
1038 update_data
1039 !> Aerosol representation id
1040 integer(kind=i_kind), intent(in) :: aero_rep_type
1041
1042 ! If an aerosol representation id has not been generated, do it now
1043 if (aero_rep_id_.eq.-1) then
1044 aero_rep_id_ = generate_int_id()
1045 end if
1046
1047 update_data%aero_rep_unique_id = aero_rep_id_
1048 update_data%aero_rep_type = int(aero_rep_type, kind=c_int)
1049 update_data%update_data = aero_rep_modal_binned_mass_create_gmd_update_data()
1050 update_data%is_malloced = .true.
1051
1052 end subroutine update_data_init_gmd
1053
1054!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1055
1056 !> Set packed update data for mode GMD
1057 subroutine update_data_set_gmd(this, section_id, GMD)
1058
1059 !> Update data
1060 class(aero_rep_update_data_modal_binned_mass_gmd_t), intent(inout) :: this
1061 !> Aerosol section id from
1062 !! camp_aero_rep_modal_binned_mass::aero_rep_modal_binned_mass_t::get_section_id
1063 integer(kind=i_kind), intent(in) :: section_id
1064 !> Updated GMD (m)
1065 real(kind=dp), intent(in) :: gmd
1066
1068 this%aero_rep_unique_id, section_id-1, gmd)
1069
1070 end subroutine update_data_set_gmd
1071
1072!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1073
1074 !> Determine the size of a binary required to pack the reaction data
1075 integer(kind=i_kind) function internal_pack_size_gmd(this, comm) &
1076 result(pack_size)
1077
1078 !> Aerosol representation update data
1079 class(aero_rep_update_data_modal_binned_mass_gmd_t), intent(in) :: this
1080 !> MPI communicator
1081 integer, intent(in) :: comm
1082
1083 pack_size = &
1084 camp_mpi_pack_size_logical(this%is_malloced, comm) + &
1085 camp_mpi_pack_size_integer(this%aero_rep_unique_id, comm)
1086
1087 end function internal_pack_size_gmd
1088
1089!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1090
1091 !> Pack the given value to the buffer, advancing position
1092 subroutine internal_bin_pack_gmd(this, buffer, pos, comm)
1093
1094 !> Aerosol representation update data
1095 class(aero_rep_update_data_modal_binned_mass_gmd_t), intent(in) :: this
1096 !> Memory buffer
1097 character, intent(inout) :: buffer(:)
1098 !> Current buffer position
1099 integer, intent(inout) :: pos
1100 !> MPI communicator
1101 integer, intent(in) :: comm
1102
1103#ifdef CAMP_USE_MPI
1104 integer :: prev_position
1105
1106 prev_position = pos
1107 call camp_mpi_pack_logical(buffer, pos, this%is_malloced, comm)
1108 call camp_mpi_pack_integer(buffer, pos, this%aero_rep_unique_id, comm)
1109 call assert(685522546, &
1110 pos - prev_position <= this%pack_size(comm))
1111#endif
1112
1113 end subroutine internal_bin_pack_gmd
1114
1115!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1116
1117 !> Unpack the given value from the buffer, advancing position
1118 subroutine internal_bin_unpack_gmd(this, buffer, pos, comm)
1119
1120 !> Aerosol representation update data
1121 class(aero_rep_update_data_modal_binned_mass_gmd_t), intent(inout) :: this
1122 !> Memory buffer
1123 character, intent(inout) :: buffer(:)
1124 !> Current buffer position
1125 integer, intent(inout) :: pos
1126 !> MPI communicator
1127 integer, intent(in) :: comm
1128
1129#ifdef CAMP_USE_MPI
1130 integer :: prev_position
1131
1132 prev_position = pos
1133 call camp_mpi_unpack_logical(buffer, pos, this%is_malloced, comm)
1134 call camp_mpi_unpack_integer(buffer, pos, this%aero_rep_unique_id, comm)
1135 call assert(855679450, &
1136 pos - prev_position <= this%pack_size(comm))
1138#endif
1139
1140 end subroutine internal_bin_unpack_gmd
1141
1142!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1143
1144 !> Finalize a GMD update data object
1145 elemental subroutine update_data_gmd_finalize(this)
1146
1147 !> Update data object to free
1148 type(aero_rep_update_data_modal_binned_mass_gmd_t), intent(inout) :: this
1149
1150 if (this%is_malloced) call aero_rep_free_update_data(this%update_data)
1151
1152 end subroutine update_data_gmd_finalize
1153
1154!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1155
1156 !> Initialize a GSD update data object
1157 subroutine update_data_init_gsd(this, update_data, aero_rep_type)
1158
1159 use camp_rand, only : generate_int_id
1160
1161 !> Aerosol representation to update
1162 class(aero_rep_modal_binned_mass_t), intent(inout) :: this
1163 !> Update data object
1164 class(aero_rep_update_data_modal_binned_mass_gsd_t), intent(out) :: &
1165 update_data
1166 !> Aerosol representation id
1167 integer(kind=i_kind), intent(in) :: aero_rep_type
1168
1169 ! If an aerosol representation id has not been generated, do it now
1170 if (aero_rep_id_.eq.-1) then
1171 aero_rep_id_ = generate_int_id()
1172 end if
1173
1174 update_data%aero_rep_unique_id = aero_rep_id_
1175 update_data%aero_rep_type = int(aero_rep_type, kind=c_int)
1176 update_data%update_data = aero_rep_modal_binned_mass_create_gsd_update_data()
1177 update_data%is_malloced = .true.
1178
1179 end subroutine update_data_init_gsd
1180
1181!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1182
1183 !> Set packed update data for mode GSD
1184 subroutine update_data_set_gsd(this, section_id, GSD)
1185
1186 !> Update data
1187 class(aero_rep_update_data_modal_binned_mass_gsd_t), intent(inout) :: this
1188 !> Aerosol section id from
1189 !! camp_aero_rep_modal_binned_mass::aero_rep_modal_binned_mass_t::get_section_id
1190 integer(kind=i_kind), intent(in) :: section_id
1191 !> Updated GSD (m)
1192 real(kind=dp), intent(in) :: gsd
1193
1195 this%aero_rep_unique_id, section_id-1, gsd)
1196
1197 end subroutine update_data_set_gsd
1198
1199!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1200
1201 !> Determine the size of a binary required to pack the reaction data
1202 integer(kind=i_kind) function internal_pack_size_gsd(this, comm) &
1203 result(pack_size)
1204
1205 !> Aerosol representation update data
1206 class(aero_rep_update_data_modal_binned_mass_gsd_t), intent(in) :: this
1207 !> MPI communicator
1208 integer, intent(in) :: comm
1209
1210 pack_size = &
1211 camp_mpi_pack_size_logical(this%is_malloced, comm) + &
1212 camp_mpi_pack_size_integer(this%aero_rep_unique_id, comm)
1213
1214 end function internal_pack_size_gsd
1215
1216!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1217
1218 !> Pack the given value to the buffer, advancing position
1219 subroutine internal_bin_pack_gsd(this, buffer, pos, comm)
1220
1221 !> Aerosol representation update data
1222 class(aero_rep_update_data_modal_binned_mass_gsd_t), intent(in) :: this
1223 !> Memory buffer
1224 character, intent(inout) :: buffer(:)
1225 !> Current buffer position
1226 integer, intent(inout) :: pos
1227 !> MPI communicator
1228 integer, intent(in) :: comm
1229
1230#ifdef CAMP_USE_MPI
1231 integer :: prev_position
1232
1233 prev_position = pos
1234 call camp_mpi_pack_logical(buffer, pos, this%is_malloced, comm)
1235 call camp_mpi_pack_integer(buffer, pos, this%aero_rep_unique_id, comm)
1236 call assert(295993259, &
1237 pos - prev_position <= this%pack_size(comm))
1238#endif
1239
1240 end subroutine internal_bin_pack_gsd
1241
1242!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1243
1244 !> Unpack the given value from the buffer, advancing position
1245 subroutine internal_bin_unpack_gsd(this, buffer, pos, comm)
1246
1247 !> Aerosol representation update data
1248 class(aero_rep_update_data_modal_binned_mass_gsd_t), intent(inout) :: this
1249 !> Memory buffer
1250 character, intent(inout) :: buffer(:)
1251 !> Current buffer position
1252 integer, intent(inout) :: pos
1253 !> MPI communicator
1254 integer, intent(in) :: comm
1255
1256#ifdef CAMP_USE_MPI
1257 integer :: prev_position
1258
1259 prev_position = pos
1260 call camp_mpi_unpack_logical(buffer, pos, this%is_malloced, comm)
1261 call camp_mpi_unpack_integer(buffer, pos, this%aero_rep_unique_id, comm)
1262 call assert(518724415, &
1263 pos - prev_position <= this%pack_size(comm))
1265#endif
1266
1267 end subroutine internal_bin_unpack_gsd
1268
1269!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1270
1271 !> Finalize a GSD update data object
1272 elemental subroutine update_data_gsd_finalize(this)
1273
1274 !> Update data object to free
1275 type(aero_rep_update_data_modal_binned_mass_gsd_t), intent(inout) :: this
1276
1277 if (this%is_malloced) call aero_rep_free_update_data(this%update_data)
1278
1279 end subroutine update_data_gsd_finalize
1280
1281!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1282
Get the size of the section of the camp_camp_state::camp_state_t::state_var array required for this a...
Initialize the aerosol representation data, validating component data and loading any required inform...
Extending-type binary pack function (Internal use only)
Extending-type binary unpack function (Internal use only)
Extending-type binary pack size (internal use only)
Get the number of Jacobian elements used in calculations of aerosol mass, volume, number,...
Get the number of instances of a specified aerosol phase.
Get the non-unique name of a chemical species by its unique name.
Get a species id on the camp_camp_state::camp_state_t::state_var array by unique name....
Get a list of unique names for each element on the camp_camp_state::camp_state_t::state_var array for...
Interface for to_string functions.
Definition util.F90:32
The abstract aero_phase_data_t structure and associated subroutines.
elemental subroutine finalize(this)
Finalize the aerosol phase data.
type(aero_phase_data_t) function, pointer constructor(phase_name, init_size)
Constructor for aero_phase_data_t.
integer(kind=i_kind) function pack_size(this, comm)
Determine the size of a binary required to pack the aerosol representation data.
The abstract aero_rep_data_t structure and associated subroutines.
The abstract aero_rep_modal_binned_mass_t structure and associated subroutines.
subroutine update_data_init_gsd(this, update_data, aero_rep_type)
Initialize a GSD update data object.
subroutine internal_bin_unpack_gsd(this, buffer, pos, comm)
Unpack the given value from the buffer, advancing position.
logical function get_section_id(this, section_name, section_id)
Get an id for a mode or bin by name for use with updates from external modules.
subroutine internal_bin_pack_gmd(this, buffer, pos, comm)
Pack the given value to the buffer, advancing position.
subroutine update_data_set_gsd(this, section_id, gsd)
Set packed update data for mode GSD.
subroutine internal_bin_unpack_gmd(this, buffer, pos, comm)
Unpack the given value from the buffer, advancing position.
subroutine internal_bin_pack_gsd(this, buffer, pos, comm)
Pack the given value to the buffer, advancing position.
integer(kind=i_kind) function internal_pack_size_gmd(this, comm)
Determine the size of a binary required to pack the reaction data.
subroutine update_data_set_gmd(this, section_id, gmd)
Set packed update data for mode GMD.
integer(kind=i_kind) function internal_pack_size_gsd(this, comm)
Determine the size of a binary required to pack the reaction data.
integer(kind=i_kind), parameter, public update_gsd
elemental subroutine update_data_gsd_finalize(this)
Finalize a GSD update data object.
integer(kind=i_kind), parameter, public update_gmd
subroutine update_data_init_gmd(this, update_data, aero_rep_type)
Initialize a GMD update object.
elemental subroutine update_data_gmd_finalize(this)
Finalize a GMD update data object.
The camp_state_t structure and associated subroutines.
Definition camp_state.F90:9
The chem_spec_data_t structure and associated subroutines.
Wrapper functions for MPI.
Definition mpi.F90:13
subroutine camp_mpi_pack_logical(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition mpi.F90:792
subroutine camp_mpi_unpack_integer(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition mpi.F90:1023
subroutine camp_mpi_unpack_logical(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition mpi.F90:1131
integer function camp_mpi_pack_size_logical(val, comm)
Determines the number of bytes required to pack the given value.
Definition mpi.F90:484
subroutine camp_mpi_pack_integer(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition mpi.F90:691
integer function camp_mpi_pack_size_integer(val, comm)
Determines the number of bytes required to pack the given value.
Definition mpi.F90:398
The property_t structure and associated subroutines.
Definition property.F90:9
Random number generators.
Definition rand.F90:9
integer(kind=i_kind) function generate_int_id()
Generate an integer id Ids will be sequential, and can only be generated by the primary process.
Definition rand.F90:435
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
Pointer type for building arrays.
Abstract aerosol representation data type.
String type for building arrays of string of various size.
Definition util.F90:38