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