CAMP 1.0.0
Chemistry Across Multiple Phases
aero_rep_single_particle.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_single_particle module.
7
8!> \page camp_aero_rep_single_particle CAMP: Single Particle Aerosol Representation
9!!
10!! The single particle aerosol representation is for use with a PartMC
11!! particle-resolved run. The \c json object for this \ref camp_aero_rep
12!! "aerosol representation" has the following format:
13!! \code{.json}
14!! { "camp-data" : [
15!! {
16!! "name" : "my single particle aero rep",
17!! "type" : "AERO_REP_SINGLE_PARTICLE"
18!! },
19!! ...
20!! ]}
21!! \endcode
22!! The key-value pair \b type is required and must be \b
23!! AERO_REP_SINGLE_PARTICLE. In this representation, particles are divided
24!! into layers. Phases in each layer are specified by the user.
25!!
26!! Phase configuration within particles follows "fractional volume overlap".
27!! The shared surface area between phases in adjacent layers is scaled by phase
28!! volume fraction within respective layers.
29!! In this configuration, the surface area of the layer interface between two phases
30!! is calculated as f_first * f_second * total_interface_surface_area where
31!! the first and second subscripts refer to the two phases in adjacent layers and
32!!
33!! f_first = volume_phase_first / volume_total_layer_first
34!! f_second = volume_phase_second / volume_total_layer_second
35!!
36!!
37!! The number concentration for each particle must be
38!! set from an external model using
39!! \c camp_aero_rep_single_particle::aero_rep_update_data_single_particle_number_t
40!! objects.
41
42!> The aero_rep_single_particle_t type and associated subroutines.
44
49 use camp_mpi
51 use camp_util, only: dp, i_kind, &
54 assert
55
56 use iso_c_binding
57
58 implicit none
59 private
60
61#define NUM_LAYERS_ this%condensed_data_int(1)
62#define AERO_REP_ID_ this%condensed_data_int(2)
63#define MAX_PARTICLES_ this%condensed_data_int(3)
64#define PARTICLE_STATE_SIZE_ this%condensed_data_int(4)
65#define NUM_INT_PROP_ 4
66#define NUM_REAL_PROP_ 0
67#define NUM_ENV_PARAM_PER_PARTICLE_ 1
68#define LAYER_PHASE_START_(l) this%condensed_data_int(NUM_INT_PROP_+l)
69#define LAYER_PHASE_END_(l) this%condensed_data_int(NUM_INT_PROP_+NUM_LAYERS_+l)
70#define TOTAL_NUM_PHASES_ (LAYER_PHASE_END_(NUM_LAYERS_))
71#define NUM_PHASES_(l) (LAYER_PHASE_END_(l)-LAYER_PHASE_START_(l)+1)
72#define PHASE_STATE_ID_(l,p) this%condensed_data_int(NUM_INT_PROP_+2*NUM_LAYERS_+LAYER_PHASE_START_(l)+p-1)
73#define PHASE_MODEL_DATA_ID_(l,p) this%condensed_data_int(NUM_INT_PROP_+2*NUM_LAYERS_+TOTAL_NUM_PHASES_+LAYER_PHASE_START_(l)+p-1)
74#define PHASE_NUM_JAC_ELEM_(l,p) this%condensed_data_int(NUM_INT_PROP_+2*NUM_LAYERS_+2*TOTAL_NUM_PHASES_+LAYER_PHASE_START_(l)+p-1)
75
76 ! Update types (These must match values in aero_rep_single_particle.c)
77 integer(kind=i_kind), parameter, public :: update_number_conc = 0
78
82
83 !> Single particle aerosol representation
84 !!
85 !! Time-invariant data related to a single particle aerosol representation.
87 !> Unique names for each instance of every chemical species in the
88 !! aerosol representaiton
89 type(string_t), allocatable, private :: unique_names_(:)
90 !> Layer names, ordered inner-most to outer-most
91 type(string_t), allocatable, private :: layer_names_(:)
92 !> First state id for the representation (only used during initialization)
93 integer(kind=i_kind) :: state_id_start = -99999
94 contains
95 !> Initialize the aerosol representation data, validating component data and
96 !! loading any required information from the \c
97 !! aero_rep_data_t::property_set. This routine should be called once for
98 !! each aerosol representation at the beginning of a model run after all
99 !! the input files have been read in. It ensures all data required during
100 !! the model run are included in the condensed data arrays.
101 procedure :: initialize
102 !> Returns the maximum number of computational particles
104 !> Initialize an update data number object
105 procedure :: update_data_initialize_number => update_data_init_number
106 !> Get the size of the section of the
107 !! \c camp_camp_state::camp_state_t::state_var array required for this
108 !! aerosol representation.
109 !!
110 !! For a single particle representation, the size will correspond to the
111 !! the sum of the sizes of a single instance of each aerosol phase
112 !! provided to \c aero_rep_single_particle::initialize()
113 procedure :: size => get_size
114 !> Get the number of state variables per-particle
115 !!
116 !! Calling functions can assume each particle has the same size on the
117 !! state array, and that individual particle states are contiguous and
118 !! arranged sequentially
119 procedure :: per_particle_size
120 !> Get a list of unique names for each element on the
121 !! \c camp_camp_state::camp_state_t::state_var array for this aerosol
122 !! representation. The list may be restricted to a particular phase and/or
123 !! aerosol species by including the phase_name and spec_name arguments.
124 !!
125 !! For a single particle representation, the unique names will be the
126 !! phase name with the species name separated by a '.'
127 procedure :: unique_names
128 !> Get a species id on the \c camp_camp_state::camp_state_t::state_var
129 !! array by its unique name. These are unique ids for each element on the
130 !! state array for this \ref camp_aero_rep "aerosol representation" and
131 !! are numbered:
132 !!
133 !! \f[x_u \in x_f ... (x_f+n-1)\f]
134 !!
135 !! where \f$x_u\f$ is the id of the element corresponding to the species
136 !! with unique name \f$u\f$ on the \c
137 !! camp_camp_state::camp_state_t::state_var array, \f$x_f\f$ is the index
138 !! of the first element for this aerosol representation on the state array
139 !! and \f$n\f$ is the total number of variables on the state array from
140 !! this aerosol representation.
141 procedure :: spec_state_id
142 !> Get the non-unique name of a species by its unique name
143 procedure :: spec_name
144 !> Get the number of instances of an aerosol phase in this representation
146 !> Get the number of Jacobian elements used in calculations of aerosol
147 !! mass, volume, number, etc. for a particular phase
148 procedure :: num_jac_elem
149 !> Returns the number of layers
150 procedure :: num_layers
151 !> Returns the number of phases in a layer or overall
152 procedure :: num_phases
153 !> Returns the number of state variables for a layer and phase
154 procedure :: phase_state_size
155 !> Returns index_pair_t type with phase_ids of adjacent phases
156 procedure :: adjacent_phases
157 !> Finalize the aerosol representation
160
161 ! Constructor for aero_rep_single_particle_t
163 procedure :: constructor
164 end interface aero_rep_single_particle_t
165
166 !> Single particle update number concentration object
167 type, extends(aero_rep_update_data_t) :: &
169 private
170 !> Flag indicating whether the update data is allocated
171 logical :: is_malloced = .false.
172 !> Unique id for finding aerosol representations during initialization
173 integer(kind=i_kind) :: aero_rep_unique_id = 0
174 !> Maximum number of computational particles
175 integer(kind=i_kind) :: maximum_computational_particles = 0
176 contains
177 !> Update the number
178 procedure :: set_number__n_m3 => update_data_set_number__n_m3
179 !> Determine the pack size of the local update data
181 !> Pack the local update data to a binary
183 !> Unpack the local update data from a binary
185 !> Finalize the number update data
188
189 !> Interface to c aerosol representation functions
190 interface
191
192 !> Allocate space for a number update
194 result(update_data) bind (c)
195 use iso_c_binding
196 !> Allocated update_data object
197 type(c_ptr) :: update_data
199
200 !> Set a new particle number concentration
202 update_data, aero_rep_unique_id, particle_id, number_conc) &
203 bind (c)
204 use iso_c_binding
205 !> Update data
206 type(c_ptr), value :: update_data
207 !> Aerosol representation unique id
208 integer(kind=c_int), value :: aero_rep_unique_id
209 !> Computational particle index
210 integer(kind=c_int), value :: particle_id
211 !> New number (m)
212 real(kind=c_double), value :: number_conc
214
215 !> Free an update data object
216 pure subroutine aero_rep_free_update_data(update_data) bind (c)
217 use iso_c_binding
218 !> Update data
219 type(c_ptr), value, intent(in) :: update_data
220 end subroutine aero_rep_free_update_data
221
222 end interface
223
224contains
225
226!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
227
228 !> Constructor for aero_rep_single_particle_t
229 function constructor() result (new_obj)
230
231 !> New aerosol representation
232 type(aero_rep_single_particle_t), pointer :: new_obj
233
234 allocate(new_obj)
235
236 end function constructor
237
238!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
239
240 !> Initialize the aerosol representation data, validating component data and
241 !! loading any required information from the \c
242 !! aero_rep_data_t::property_set. This routine should be called once for
243 !! each aerosol representation at the beginning of a model run after all
244 !! the input files have been read in. It ensures all data required during
245 !! the model run are included in the condensed data arrays.
246 subroutine initialize(this, aero_phase_set, spec_state_id)
247
248 !> Aerosol representation data
249 class(aero_rep_single_particle_t), intent(inout) :: this
250 !> The set of aerosol phases
251 type(aero_phase_data_ptr), pointer, intent(in) :: aero_phase_set(:)
252 !> Beginning state id for this aerosol representation in the model species
253 !! state array
254 integer(kind=i_kind), intent(in) :: spec_state_id
255
256 ! Unordered layer names (only used during initialization)
257 type(string_t), allocatable :: layer_names_unordered(:)
258 ! Cover names (only used during initialization)
259 type(string_t), allocatable :: cover_names_unordered(:)
260 ! Index in layer_names_unordered for each layer from inner- to outer-most
261 ! layer
262 integer(kind=i_kind), allocatable :: ordered_layer_id(:)
263
264 character(len=:), allocatable :: key_name, layer_name, layer_covers, &
265 phase_name
266 type(property_t), pointer :: layers, layer
267 type(property_ptr), allocatable :: phases(:)
268 integer(kind=i_kind) :: i_particle, i_phase, i_layer, i_aero, curr_id
269 integer(kind=i_kind) :: i_cover, j_phase, j_layer, i_map, curr_phase
270 integer(kind=i_kind) :: num_phases, num_int_param, num_float_param, &
271 num_particles
272 logical :: found
273
274 ! Get the maximum number of computational particles
275 key_name = "maximum computational particles"
276 call assert_msg(331697425, &
277 this%property_set%get_int(key_name, num_particles), &
278 "Missing maximum number of computational particles")
279
280 ! Get the set of layers
281 key_name = "layers"
282 call assert_msg(314292954, &
283 this%property_set%get_property_t(key_name, layers), &
284 "Missing layers for single-particle aerosol "// &
285 "representation '"//this%rep_name//"'")
286 call assert_msg(168669831, layers%size() .gt. 0, &
287 "No Layers specified for single-particle layer "// &
288 "aerosol representation '"//this%rep_name//"'")
289
290 ! Allocate space for the working arrays
291 allocate(phases(layers%size()))
292 allocate(cover_names_unordered(layers%size()))
293 allocate(layer_names_unordered(layers%size()))
294
295 ! Loop through the layers, adding names and counting the spaces needed
296 ! on the condensed data arrays, and counting the total phases instances
297 num_phases = 0
298 call layers%iter_reset()
299 do i_layer = 1, layers%size()
300
301 ! Get the layer properties
302 call assert_msg(303808978, layers%get_property_t(val=layer), &
303 "Invalid structure for layer '"// &
304 layer_names_unordered(i_layer)%string// &
305 "' in single-particle layer representation '"// &
306 this%rep_name//"'")
307
308 ! Get the layer name
309 key_name = "name"
310 call assert_msg(364496472, layer%get_string(key_name, layer_name), &
311 "Missing layer name in single-particle layer aerosol "// &
312 "representation '"//this%rep_name//"'")
313 layer_names_unordered(i_layer)%string = layer_name
314
315 ! Get the cover name
316 key_name = "covers"
317 call assert_msg(350939595, layer%get_string(key_name, layer_covers), &
318 "Missing cover name in layer'"// &
319 layer_names_unordered(i_layer)%string// &
320 "' in single-particle layer aerosol representation '"// &
321 this%rep_name//"'")
322 cover_names_unordered(i_layer)%string = layer_covers
323
324 ! Get the set of phases
325 key_name = "phases"
326 call assert_msg(647756433, &
327 layer%get_property_t(key_name, phases(i_layer)%val_), &
328 "Missing phases for layer '"// &
329 layer_names_unordered(i_layer)%string// &
330 "' in single-particle layer aerosol representation '"// &
331 this%rep_name//"'")
332
333 ! Add the phases to the counter
334 call assert_msg(002679882, phases(i_layer)%val_%size().gt.0, &
335 "No phases specified for layer '"// &
336 layer_names_unordered(i_layer)%string// &
337 "' in single-particle layer aerosol representation '"// &
338 this%rep_name//"'")
339
340 ! add to running total of phase count
341 num_phases = num_phases + phases(i_layer)%val_%size()
342
343 call layers%iter_next()
344 end do
345
346 ! get the map of layer names after reordering
347 ordered_layer_id = ordered_layer_ids(layer_names_unordered, &
348 cover_names_unordered)
349
350 ! set the layer names
351 allocate(this%layer_names_(size(ordered_layer_id)))
352 this%layer_names_(:) = layer_names_unordered(ordered_layer_id(:))
353
354 ! Allocate condensed data arrays
355 num_int_param = num_int_prop_ + 2 * layers%size() + 3 * num_phases
356 num_float_param = num_real_prop_
357 allocate(this%condensed_data_int(num_int_param))
358 allocate(this%condensed_data_real(num_float_param))
359 this%condensed_data_int(:) = int(0, kind=i_kind)
360 this%condensed_data_real(:) = real(0.0, kind=dp)
361
362 ! Save space for the environment-dependent parameters
363 this%num_env_params = num_env_param_per_particle_ * num_particles
364
365 ! Save representation dimensions
366 num_layers_ = layers%size()
367 max_particles_ = num_particles
368
369 ! validate phase names, assign aero_phase pointers for each phase in
370 ! each layer in each particle, and set PHASE_STATE_ID and
371 ! PHASE_MODEL_DATA_ID for each phase
372 allocate(this%aero_phase(num_phases * num_particles))
373 allocate(this%aero_phase_is_at_surface(num_phases * num_particles))
374 curr_phase = 1
375 do i_layer = 1, size(ordered_layer_id)
376 j_layer = ordered_layer_id(i_layer)
377
378 ! Set the starting and ending indices for the phases in this layer
379 layer_phase_start_(i_layer) = curr_phase
380 layer_phase_end_(i_layer) = curr_phase + phases(j_layer)%val_%size() - 1
381
382 curr_phase = curr_phase + phases(j_layer)%val_%size()
383 end do
384
385 curr_id = spec_state_id
386 this%state_id_start = spec_state_id
387 curr_phase = 1
388 do i_layer = 1, size(ordered_layer_id)
389 j_layer = ordered_layer_id(i_layer)
390 ! Loop through the phases and make sure they exist
391 call phases(j_layer)%val_%iter_reset()
392 do i_phase = 1, phases(j_layer)%val_%size()
393
394 ! Get the phase name
395 call assert_msg(566480284, &
396 phases(j_layer)%val_%get_string(val=phase_name), &
397 "Non-string phase name for layer '"// &
398 layer_names_unordered(j_layer)%string// &
399 "' in single-particle layer aerosol representation '"// &
400 this%rep_name//"'")
401 ! find phase and set pointer and indices
402 found = .false.
403 do j_phase = 1, size(aero_phase_set)
404 if (aero_phase_set(j_phase)%val%name() .eq. phase_name) then
405 found = .true.
406 do i_particle = 0, num_particles-1
407 this%aero_phase(i_particle*num_phases + curr_phase) = &
408 aero_phase_set(j_phase)
409 if (i_layer .eq. num_layers_) then
410 this%aero_phase_is_at_surface(i_particle*num_phases + curr_phase) = &
411 .true.
412 else
413 this%aero_phase_is_at_surface(i_particle*num_phases + curr_phase) = &
414 .false.
415 end if
416 end do
417 phase_state_id_(i_layer,i_phase) = curr_id
418 phase_model_data_id_(i_layer,i_phase) = j_phase
419 curr_id = curr_id + aero_phase_set(j_phase)%val%size()
420 curr_phase = curr_phase + 1
421 exit
422 end if
423 end do
424 call assert(373124707, found)
425
426 call phases(j_layer)%val_%iter_next()
427 end do
428 end do
429 particle_state_size_ = curr_id - spec_state_id
430
431 ! Initialize the aerosol representation id
432 aero_rep_id_ = -1
433
434 ! Set the unique names for the chemical species
435 this%unique_names_ = this%unique_names( )
436
437 end subroutine initialize
438
439!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
440
441 !> Returns the maximum nunmber of computational particles
442 integer(kind=i_kind) function maximum_computational_particles(this)
443
444 !> Aerosol representation data
445 class(aero_rep_single_particle_t), intent(in) :: this
446
447 maximum_computational_particles = max_particles_
448
450
451!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
452
453 !> Get the size of the section of the
454 !! \c camp_camp_state::camp_state_t::state_var array required for this
455 !! aerosol representation.
456 !!
457 !! For a single particle representation, the size will correspond to the
458 !! the sum of the sizes of a single instance of each aerosol phase
459 !! provided to \c aero_rep_single_particle::initialize()
460 function get_size(this) result (state_size)
461
462 !> Size on the state array
463 integer(kind=i_kind) :: state_size
464 !> Aerosol representation data
465 class(aero_rep_single_particle_t), intent(in) :: this
466
467 state_size = max_particles_ * particle_state_size_
468
469 end function get_size
470
471!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
472
473 !> Get the number of state variables per-particle
474 !!
475 !! Calling functions can assume each particle has the same size on the
476 !! state array, and that individual particle states are contiguous and
477 !! arranged sequentially
478 function per_particle_size(this) result(state_size)
479
480 !> Size on the state array per particle
481 integer(kind=i_kind) :: state_size
482 !> Aerosol representation data
483 class(aero_rep_single_particle_t), intent(in) :: this
484
485 state_size = particle_state_size_
486
487 end function per_particle_size
488
489!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
490
491 !> Get a list of unique names for each element on the
492 !! \c camp_camp_state::camp_state_t::state_var array for this aerosol
493 !! representation. The list may be restricted to a particular phase and/or
494 !! aerosol species by including the phase_name and spec_name arguments.
495 !!
496 !! For a single particle representation, the unique names will be a 'P'
497 !! followed by the computational particle number, a '.', the phase name,
498 !! another '.', and the species name.
499 function unique_names(this, phase_name, tracer_type, spec_name, &
500 phase_is_at_surface)
501
502 use camp_util, only : integer_to_string
503 !> List of unique names
504 type(string_t), allocatable :: unique_names(:)
505 !> Aerosol representation data
506 class(aero_rep_single_particle_t), intent(in) :: this
507 !> Aerosol phase name
508 character(len=*), optional, intent(in) :: phase_name
509 !> Aerosol-phase species tracer type
510 integer(kind=i_kind), optional, intent(in) :: tracer_type
511 !> Aerosol-phase species name
512 character(len=*), optional, intent(in) :: spec_name
513 !> Indicates if aerosol phase is at the surface of particle
514 logical, optional, intent(in) :: phase_is_at_surface
515
516 integer :: i_particle, i_layer, i_phase, i_spec, j_spec
517 integer :: num_spec, curr_tracer_type
518 type(string_t), allocatable :: spec_names(:)
519
520 ! copy saved unique names when available and no filters are included
521 if (.not. present(phase_name) .and. &
522 .not. present(tracer_type) .and. &
523 .not. present(spec_name) .and. &
524 allocated(this%unique_names_)) then
525 unique_names = this%unique_names_
526 return
527 end if
528
529 ! count the number of unique names
530 num_spec = 0
531 do i_phase = 1, size(this%aero_phase)
532 if (present(phase_name)) then
533 if(phase_name .ne. this%aero_phase(i_phase)%val%name()) cycle
534 end if
535 if (present(phase_is_at_surface)) then
536 if (phase_is_at_surface .neqv. &
537 this%aero_phase_is_at_surface(i_phase)) cycle
538 end if
539 if (present(spec_name) .or. present(tracer_type)) then
540 spec_names = this%aero_phase(i_phase)%val%get_species_names()
541 do j_spec = 1, size(spec_names)
542 if (present(spec_name)) then
543 if (spec_name .ne. spec_names(j_spec)%string) cycle
544 end if
545 if (present(tracer_type)) then
546 curr_tracer_type = &
547 this%aero_phase(i_phase)%val%get_species_type( &
548 spec_names(j_spec)%string)
549 if (tracer_type .ne. curr_tracer_type) cycle
550 end if
551 num_spec = num_spec + 1
552 end do
553 deallocate(spec_names)
554 else
555 num_spec = num_spec + this%aero_phase(i_phase)%val%size()
556 end if
557 end do
558
559 ! allocate space for the unique names and assign them
560 num_spec = num_spec / max_particles_ ! we need per-particle value for indexing
561 allocate(unique_names(num_spec*max_particles_))
562 i_spec = 1
563 do i_layer = 1, num_layers_
564 do i_phase = layer_phase_start_(i_layer), layer_phase_end_(i_layer)
565 if (present(phase_name)) then
566 if(phase_name .ne. this%aero_phase(i_phase)%val%name()) cycle
567 end if
568 if (present(phase_is_at_surface)) then
569 if (phase_is_at_surface .neqv. &
570 this%aero_phase_is_at_surface(i_phase)) cycle
571 end if
572 spec_names = this%aero_phase(i_phase)%val%get_species_names()
573 do j_spec = 1, this%aero_phase(i_phase)%val%size()
574 if (present(spec_name)) then
575 if (spec_name .ne. spec_names(j_spec)%string) cycle
576 end if
577 if (present(tracer_type)) then
578 curr_tracer_type = &
579 this%aero_phase(i_phase)%val%get_species_type( &
580 spec_names(j_spec)%string)
581 if (tracer_type .ne. curr_tracer_type) cycle
582 end if
583 do i_particle = 1, max_particles_
584 unique_names((i_particle-1)*num_spec+i_spec)%string = 'P'// &
585 trim(integer_to_string(i_particle))//"."// &
586 this%layer_names_(i_layer)%string//"."// &
587 this%aero_phase(i_phase)%val%name()//"."// &
588 spec_names(j_spec)%string
589 end do
590 i_spec = i_spec + 1
591 end do
592 deallocate(spec_names)
593 end do
594 end do
595
596 end function unique_names
597
598!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
599
600 !> Get a species id on the \c camp_camp_state::camp_state_t::state_var
601 !! array by its unique name. These are unique ids for each element on the
602 !! state array for this \ref camp_aero_rep "aerosol representation" and
603 !! are numbered:
604 !!
605 !! \f[x_u \in x_f ... (x_f+n-1)\f]
606 !!
607 !! where \f$x_u\f$ is the id of the element corresponding to the species
608 !! with unique name \f$u\f$ on the \c
609 !! camp_camp_state::camp_state_t::state_var array, \f$x_f\f$ is the index
610 !! of the first element for this aerosol representation on the state array
611 !! and \f$n\f$ is the total number of variables on the state array from
612 !! this aerosol representation.
613 function spec_state_id(this, unique_name) result (spec_id)
614
615 !> Species state id
616 integer(kind=i_kind) :: spec_id
617 !> Aerosol representation data
618 class(aero_rep_single_particle_t), intent(in) :: this
619 !> Unique name
620 character(len=*), intent(in) :: unique_name
621
622 integer(kind=i_kind) :: i_spec
623
624 spec_id = this%state_id_start
625 do i_spec = 1, size(this%unique_names_)
626 if (this%unique_names_(i_spec)%string .eq. unique_name) then
627 return
628 end if
629 spec_id = spec_id + 1
630 end do
631 call die_msg( 449087541, "Cannot find species '"//unique_name//"'" )
632
633 end function spec_state_id
634
635!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
636
637 !> Get the non-unique name of a species in this aerosol representation by
638 !! id.
639 function spec_name(this, unique_name)
640
641 !> Chemical species name
642 character(len=:), allocatable :: spec_name
643 !> Aerosol representation data
644 class(aero_rep_single_particle_t), intent(in) :: this
645 !> Unique name of the species in this aerosol representation
646 character(len=*), intent(in) :: unique_name
647
648 type(string_t) :: l_unique_name
649 type(string_t), allocatable :: substrs(:)
650
651 l_unique_name%string = unique_name
652 substrs = l_unique_name%split(".")
653 call assert(407537518, size( substrs ) .eq. 4 )
654 spec_name = substrs(4)%string
655
656 end function spec_name
657
658!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
659
660 !> Get the number of instances of a specified aerosol phase. In the single
661 !! particle representation with layers, a phase can exist in multiple layers
662 !! in one particle.
663 integer(kind=i_kind) function num_phase_instances(this, phase_name, &
664 is_at_surface)
665
666 !> Aerosol representation data
667 class(aero_rep_single_particle_t), intent(in) :: this
668 !> Aerosol phase name
669 character(len=*), intent(in) :: phase_name
670 !> Indicates if aerosol phase is at the surface of particle
671 logical, intent(in), optional :: is_at_surface
672
673 integer(kind=i_kind) :: i_phase, i_layer, phase_index
674
676 if (present(is_at_surface)) then
677 do i_layer = 1, num_layers_
678 do i_phase = layer_phase_start_(i_layer), layer_phase_end_(i_layer)
679 if (this%aero_phase(i_phase)%val%name() .eq. phase_name) then
680 if (this%aero_phase_is_at_surface(i_phase) .eqv. is_at_surface) then
682 end if
683 end if
684 end do
685 end do
686 else
687 do i_layer = 1, num_layers_
688 do i_phase = layer_phase_start_(i_layer), layer_phase_end_(i_layer)
689 if (this%aero_phase(i_phase)%val%name() .eq. phase_name) then
691 end if
692 end do
693 end do
694 end if
696
697 end function num_phase_instances
698
699!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
700
701 !> Get the number of Jacobian elements used in calculations of aerosol mass,
702 !! volume, number, etc. for a particular phase
703 function num_jac_elem(this, phase_id)
704
705 !> Number of Jacobian elements used
706 integer(kind=i_kind) :: num_jac_elem
707 !> Aerosol respresentation data
708 class(aero_rep_single_particle_t), intent(in) :: this
709 !> Aerosol phase id
710 integer(kind=i_kind), intent(in) :: phase_id
711
712 integer(kind=i_kind) :: i_phase
713
714 call assert_msg(927040495, phase_id .ge. 1 .and. &
715 phase_id .le. size( this%aero_phase ), &
716 "Aerosol phase index out of range. Got "// &
717 trim( integer_to_string( phase_id ) )//", expected 1:"// &
718 trim( integer_to_string( size( this%aero_phase ) ) ) )
719 num_jac_elem = 0
720 do i_phase = 1, total_num_phases_
722 this%aero_phase(i_phase)%val%num_jac_elem()
723 end do
724
725 end function num_jac_elem
726
727!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
728
729 !> Returns the number of layers
730 integer function num_layers(this)
731
732 !> Aerosol representation data
733 class(aero_rep_single_particle_t), intent(in) :: this
734
735 num_layers = num_layers_
736
737 end function num_layers
738
739!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
740
741 !> Returns the number of phases in a layer or overall
742 integer function num_phases(this, layer)
743
744 !> Aerosol representation data
745 class(aero_rep_single_particle_t), intent(in) :: this
746 !> Layer id
747 integer, optional, intent(in) :: layer
748
749 if (present(layer)) then
750 num_phases = layer_phase_end_(layer) - layer_phase_start_(layer) + 1
751 else
752 num_phases = layer_phase_end_(num_layers_) - layer_phase_start_(1) + 1
753 end if
754
755 end function num_phases
756
757!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
758
759 !> Returns the number of state variables for a layer and phase
760 integer function phase_state_size(this, layer, phase)
761
762 use camp_util, only : die_msg
763
764 !> Aerosol representation data
765 class(aero_rep_single_particle_t), intent(in) :: this
766 !> Layer id
767 integer, optional, intent(in) :: layer
768 !> Phase id
769 integer, optional, intent(in) :: phase
770
771 if (present(layer) .and. present(phase)) then
772 if (layer .eq. num_layers_ .and. phase .eq. num_phases_(layer)) then
773 phase_state_size = phase_state_id_(1,1) + particle_state_size_ - &
774 phase_state_id_(layer, phase)
775 else if (phase .eq. num_phases_(layer)) then
776 phase_state_size = phase_state_id_(layer+1, 1) - &
777 phase_state_id_(layer, phase)
778 else
779 phase_state_size = phase_state_id_(layer, phase+1) - &
780 phase_state_id_(layer, phase)
781 end if
782 else if (present(layer)) then
783 if (layer .eq. num_layers_) then
784 phase_state_size = phase_state_id_(1,1) + particle_state_size_ - &
785 phase_state_id_(layer, 1)
786 else
787 phase_state_size = phase_state_id_(layer+1, 1) - &
788 phase_state_id_(layer, 1)
789 end if
790 else if (present(phase)) then
791 call die_msg(917793122, "Must specify layer if including phase is "// &
792 "state size request")
793 else
794 phase_state_size = particle_state_size_
795 end if
796
797 end function phase_state_size
798
799!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
800
801 !> Determine is specified phase(s) exist in adjacent layers. Returns array
802 !! of phase_ids for adjacent phases first and second.
803
804 function adjacent_phases(this, phase_name_first, &
805 phase_name_second) result(index_pairs)
806
807 !> Aerosol representation data
808 class(aero_rep_single_particle_t), intent(in) :: this
809 character(len=*), intent(in) :: phase_name_first
810 character(len=*), intent(in) :: phase_name_second
811 type(index_pair_t), allocatable :: temp_index_pairs(:), index_pairs(:)
812
813 integer(kind=i_kind), allocatable :: layer_first(:)
814 integer(kind=i_kind), allocatable :: layer_second(:)
815 integer(kind=i_kind), allocatable :: phase_id_first_(:)
816 integer(kind=i_kind), allocatable :: phase_id_second_(:)
817 integer(kind=i_kind) :: i_layer, i_phase, i_instance, j_phase
818
819 allocate(layer_first(num_layers_))
820 allocate(layer_second(num_layers_))
821 allocate(phase_id_first_(num_layers_))
822 allocate(phase_id_second_(num_layers_))
823 layer_first = -9999
824 layer_second = -9999
825 phase_id_first_ = -9999
826 phase_id_second_ = -9999
827
828 ! Loop over layers and phases to determine where each input phase exists
829 i_instance = 1
830 do i_layer = 1, num_layers_
831 do i_phase = 1, num_phases_(i_layer)
832 if (phase_name_first .eq. phase_name_second) then
833 if (this%aero_phase(i_instance)%val%name() .eq. phase_name_first) then
834 layer_first(i_layer) = phase_model_data_id_(i_layer, i_phase)
835 phase_id_first_(i_layer) = i_instance
836 end if
837 else
838 if (this%aero_phase(i_instance)%val%name() .eq. phase_name_first) then
839 layer_first(i_layer) = phase_model_data_id_(i_layer, i_phase)
840 phase_id_first_(i_layer) = i_instance
841 else if (this%aero_phase(i_instance)%val%name() .eq. phase_name_second) then
842 layer_second(i_layer) = phase_model_data_id_(i_layer, i_phase)
843 phase_id_second_(i_layer) = i_instance
844 end if
845 end if
846 i_instance = i_instance + 1
847 end do
848 end do
849
850 ! Find out where the pairs exist and assign to temp_index_pairs
851 allocate(temp_index_pairs(i_instance))
852 i_instance = 1
853 do i_layer = 1, num_layers_-1
854 do i_phase = 1, num_phases_(i_layer)
855 do j_phase = 1, num_phases_(i_layer+1)
856 if (phase_name_first .eq. phase_name_second) then
857 if (layer_first(i_layer) .eq. phase_model_data_id_(i_layer,i_phase) .and. &
858 layer_first(i_layer+1) .eq. phase_model_data_id_(i_layer+1,j_phase)) then
859 temp_index_pairs(i_instance)%first_ = phase_id_first_(i_layer)
860 temp_index_pairs(i_instance)%second_ = phase_id_first_(i_layer+1)
861 i_instance = i_instance + 1
862 end if
863 else
864 if (layer_first(i_layer) .eq. phase_model_data_id_(i_layer, i_phase) .and. &
865 layer_second(i_layer+1) .eq. phase_model_data_id_(i_layer+1, j_phase)) then
866 temp_index_pairs(i_instance)%first_ = phase_id_first_(i_layer)
867 temp_index_pairs(i_instance)%second_ = phase_id_second_(i_layer+1)
868 i_instance = i_instance + 1
869 else if (layer_second(i_layer) .eq. phase_model_data_id_(i_layer, i_phase) .and. &
870 layer_first(i_layer+1) .eq. phase_model_data_id_(i_layer+1, j_phase)) then
871 temp_index_pairs(i_instance)%first_ = phase_id_second_(i_layer)
872 temp_index_pairs(i_instance)%second_ = phase_id_first_(i_layer+1)
873 i_instance = i_instance + 1
874 end if
875 end if
876 end do
877 end do
878 end do
879
880 allocate(index_pairs(i_instance-1))
881 index_pairs(:)%first_ = temp_index_pairs(1:i_instance-1)%first_
882 index_pairs(:)%second_ = temp_index_pairs(1:i_instance-1)%second_
883 deallocate(temp_index_pairs)
884
885 end function adjacent_phases
886
887!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
888
889 !> Finalize the aerosol representation
890 subroutine finalize(this)
891
892 !> Aerosol representation data
893 type(aero_rep_single_particle_t), intent(inout) :: this
894
895 if (allocated(this%rep_name)) deallocate(this%rep_name)
896 if (allocated(this%aero_phase)) then
897 ! The core will deallocate the aerosol phases
898 call this%aero_phase(:)%dereference()
899 deallocate(this%aero_phase)
900 end if
901 if (allocated(this%unique_names_)) deallocate(this%unique_names_)
902 if (allocated(this%layer_names_)) deallocate(this%layer_names_)
903 if (associated(this%property_set)) deallocate(this%property_set)
904 if (allocated(this%condensed_data_real)) &
905 deallocate(this%condensed_data_real)
906 if (allocated(this%condensed_data_int)) &
907 deallocate(this%condensed_data_int)
908
909 end subroutine finalize
910
911!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
912
913 !> Finalize the aerosol representation array
914 subroutine finalize_array(this_array)
915
916 !> Aerosol representation data
917 type(aero_rep_single_particle_t), intent(inout) :: this_array(:)
918
919 integer(kind=i_kind) :: i_aero
920
921 do i_aero = 1, size(this_array)
922 call finalize(this_array(i_aero))
923 end do
924
925 end subroutine finalize_array
926
927!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
928
929 !> Initialize an update data object
930 subroutine update_data_init_number(this, update_data, aero_rep_type)
931
932 use camp_rand, only : generate_int_id
933
934 !> Aerosol representation to update
935 class(aero_rep_single_particle_t), intent(inout) :: this
936 !> Update data object
937 class(aero_rep_update_data_single_particle_number_t), intent(out) :: &
938 update_data
939 !> Aerosol representaiton id
940 integer(kind=i_kind), intent(in) :: aero_rep_type
941
942 ! If an aerosol representation id has not been generated, do it now
943 if (aero_rep_id_.eq.-1) then
944 aero_rep_id_ = generate_int_id()
945 end if
946
947 update_data%aero_rep_unique_id = aero_rep_id_
948 update_data%maximum_computational_particles = &
949 this%maximum_computational_particles( )
950 update_data%aero_rep_type = int(aero_rep_type, kind=c_int)
951 update_data%update_data = &
953 update_data%is_malloced = .true.
954
955 end subroutine update_data_init_number
956
957!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
958
959 !> Order layer array from inner most layer to outermost
960 function ordered_layer_ids(layer_names_unordered, cover_names_unordered)
961
962 !> Layer names in original order
963 type(string_t), intent(in) :: layer_names_unordered(:)
964 !> Name of "covered" layer for each layer in layer_name_unordered
965 type(string_t), intent(in) :: cover_names_unordered(:)
966 !> Index of name in layer_name_unordered for each layer after reordering
967 integer, allocatable :: ordered_layer_ids(:)
968
969 integer(kind=i_kind) :: i_layer, j_layer, i_cover
970
971 ! Ensure layer names do not repeat
972 do i_layer = 1, size(layer_names_unordered)
973 do j_layer = 1, size(layer_names_unordered)
974 if (i_layer .eq. j_layer) cycle
975 call assert_msg(781626922, layer_names_unordered(i_layer)%string .ne. &
976 layer_names_unordered(j_layer)%string, &
977 "Duplicate layer name in single particle "// &
978 "representation: '"// &
979 trim(layer_names_unordered(i_layer)%string)//"'")
980 end do
981 end do
982
983 allocate(ordered_layer_ids(size(layer_names_unordered)))
984
985 ! Search for innermost layer with cover set to 'none'
986 do i_layer = 1, size(layer_names_unordered)
987 if (cover_names_unordered(i_layer)%string == "none") then
988 ordered_layer_ids(1) = i_layer
989 end if
990 end do
991
992 ! Assign each layer working outwards from center of particle
993 do i_cover = 2, size(ordered_layer_ids)
994 do i_layer = 1, size(layer_names_unordered)
995 if (layer_names_unordered(ordered_layer_ids(i_cover-1))%string &
996 .eq. cover_names_unordered(i_layer)%string) then
997 ordered_layer_ids(i_cover) = i_layer
998 exit
999 end if
1000 end do
1001 end do
1002
1003 end function ordered_layer_ids
1004
1005!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1006
1007 !> Set packed update data for particle number (#/m3) for a particular
1008 !! computational particle.
1009 subroutine update_data_set_number__n_m3(this, particle_id, number_conc)
1010
1011 !> Update data
1012 class(aero_rep_update_data_single_particle_number_t), intent(inout) :: &
1013 this
1014 !> Computational particle index
1015 integer(kind=i_kind), intent(in) :: particle_id
1016 !> Updated number
1017 real(kind=dp), intent(in) :: number_conc
1018
1019 call assert_msg(611967802, this%is_malloced, &
1020 "Trying to set number of uninitialized update object.")
1021 call assert_msg(689085496, particle_id .ge. 1 .and. &
1022 particle_id .le. this%maximum_computational_particles, &
1023 "Invalid computational particle index: "// &
1024 trim(integer_to_string(particle_id)))
1026 this%get_data(), this%aero_rep_unique_id, particle_id-1, &
1027 number_conc)
1028
1029 end subroutine update_data_set_number__n_m3
1030
1031!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1032
1033 !> Determine the size of a binary required to pack the reaction data
1034 integer(kind=i_kind) function internal_pack_size_number(this, comm) &
1035 result(pack_size)
1036
1037 !> Aerosol representation update data
1038 class(aero_rep_update_data_single_particle_number_t), intent(in) :: this
1039 !> MPI communicator
1040 integer, intent(in) :: comm
1041
1042 pack_size = &
1043 camp_mpi_pack_size_logical(this%is_malloced, comm) + &
1044 camp_mpi_pack_size_integer(this%maximum_computational_particles, comm) + &
1045 camp_mpi_pack_size_integer(this%aero_rep_unique_id, comm)
1046
1047 end function internal_pack_size_number
1048
1049!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1050
1051 !> Pack the given value to the buffer, advancing position
1052 subroutine internal_bin_pack_number(this, buffer, pos, comm)
1053
1054 !> Aerosol representation update data
1055 class(aero_rep_update_data_single_particle_number_t), intent(in) :: this
1056 !> Memory buffer
1057 character, intent(inout) :: buffer(:)
1058 !> Current buffer position
1059 integer, intent(inout) :: pos
1060 !> MPI communicator
1061 integer, intent(in) :: comm
1062
1063#ifdef CAMP_USE_MPI
1064 integer :: prev_position
1065
1066 prev_position = pos
1067 call camp_mpi_pack_logical(buffer, pos, this%is_malloced, comm)
1068 call camp_mpi_pack_integer(buffer, pos, &
1069 this%maximum_computational_particles, comm)
1070 call camp_mpi_pack_integer(buffer, pos, this%aero_rep_unique_id, comm)
1071 call assert(411585487, &
1072 pos - prev_position <= this%pack_size(comm))
1073#endif
1074
1075 end subroutine internal_bin_pack_number
1076
1077!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1078
1079 !> Unpack the given value from the buffer, advancing position
1080 subroutine internal_bin_unpack_number(this, buffer, pos, comm)
1081
1082 !> Aerosol representation update data
1083 class(aero_rep_update_data_single_particle_number_t), intent(inout) :: this
1084 !> Memory buffer
1085 character, intent(inout) :: buffer(:)
1086 !> Current buffer position
1087 integer, intent(inout) :: pos
1088 !> MPI communicator
1089 integer, intent(in) :: comm
1090
1091#ifdef CAMP_USE_MPI
1092 integer :: prev_position
1093
1094 prev_position = pos
1095 call camp_mpi_unpack_logical(buffer, pos, this%is_malloced, comm)
1096 call camp_mpi_unpack_integer(buffer, pos, &
1097 this%maximum_computational_particles, comm)
1098 call camp_mpi_unpack_integer(buffer, pos, this%aero_rep_unique_id, comm)
1099 call assert(351557153, &
1100 pos - prev_position <= this%pack_size(comm))
1102#endif
1103
1104 end subroutine internal_bin_unpack_number
1105
1106!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1107
1108 !> Finalize a number update data object
1110
1111 !> Update data object to free
1112 type(aero_rep_update_data_single_particle_number_t), intent(inout) :: this
1113
1114 if (this%is_malloced) call aero_rep_free_update_data(this%update_data)
1115
1116 end subroutine update_data_number_finalize
1117
1118!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1119
1120 !> Finalize an array of number update data objects
1122
1123 !> Update data objects to free
1124 type(aero_rep_update_data_single_particle_number_t), intent(inout) :: &
1125 this(:)
1126
1127 integer(kind=i_kind) :: i_aero
1128
1129 do i_aero = 1, size(this)
1130 call update_data_number_finalize(this(i_aero))
1131 end do
1132
1134
1135!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1136
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.
subroutine finalize_array(this)
Finalize the aerosol phase data.
type(aero_phase_data_t) function, pointer constructor(phase_name, init_size)
Constructor for aero_phase_data_t.
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.
The abstract aero_rep_data_t structure and associated subroutines.
The aero_rep_single_particle_t type and associated subroutines.
integer(kind=i_kind) function internal_pack_size_number(this, comm)
Determine the size of a binary required to pack the reaction data.
subroutine update_data_init_number(this, update_data, aero_rep_type)
Initialize an update data object.
integer function num_layers(this)
Returns the number of layers.
subroutine update_data_number_finalize_array(this)
Finalize an array of number update data objects.
subroutine internal_bin_unpack_number(this, buffer, pos, comm)
Unpack the given value from the buffer, advancing position.
integer(kind=i_kind) function per_particle_size(this)
Get the number of state variables per-particle.
integer(kind=i_kind), parameter, public update_number_conc
subroutine update_data_number_finalize(this)
Finalize a number update data object.
subroutine update_data_set_number__n_m3(this, particle_id, number_conc)
Set packed update data for particle number (#/m3) for a particular computational particle.
integer function, dimension(:), allocatable, public ordered_layer_ids(layer_names_unordered, cover_names_unordered)
Order layer array from inner most layer to outermost.
subroutine internal_bin_pack_number(this, buffer, pos, comm)
Pack the given value to the buffer, advancing position.
type(index_pair_t) function, dimension(:), allocatable adjacent_phases(this, phase_name_first, phase_name_second)
Determine is specified phase(s) exist in adjacent layers. Returns array of phase_ids for adjacent pha...
integer function phase_state_size(this, layer, phase)
Returns the number of state variables for a layer and phase.
integer function num_phases(this, layer)
Returns the number of phases in a layer or overall.
integer(kind=i_kind) function maximum_computational_particles(this)
Returns the maximum nunmber of computational particles.
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
character(len=camp_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition util.F90:839
Pointer type for building arrays.
Abstract aerosol representation data type.
Define index_pair array for adjacent_phases functions.
String type for building arrays of string of various size.
Definition util.F90:53