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 !> Finalize the aerosol representation
146 final :: finalize
148
149 ! Constructor for aero_rep_single_particle_t
151 procedure :: constructor
152 end interface aero_rep_single_particle_t
153
154 !> Single particle update number concentration object
155 type, extends(aero_rep_update_data_t) :: &
157 private
158 !> Flag indicating whether the update data is allocated
159 logical :: is_malloced = .false.
160 !> Unique id for finding aerosol representations during initialization
161 integer(kind=i_kind) :: aero_rep_unique_id = 0
162 !> Maximum number of computational particles
163 integer(kind=i_kind) :: maximum_computational_particles = 0
164 contains
165 !> Update the number
166 procedure :: set_number__n_m3 => update_data_set_number__n_m3
167 !> Determine the pack size of the local update data
169 !> Pack the local update data to a binary
171 !> Unpack the local update data from a binary
173 !> Finalize the number update data
176
177 !> Interface to c aerosol representation functions
178 interface
179
180 !> Allocate space for a number update
182 result(update_data) bind (c)
183 use iso_c_binding
184 !> Allocated update_data object
185 type(c_ptr) :: update_data
187
188 !> Set a new particle number concentration
190 update_data, aero_rep_unique_id, particle_id, number_conc) &
191 bind (c)
192 use iso_c_binding
193 !> Update data
194 type(c_ptr), value :: update_data
195 !> Aerosol representation unique id
196 integer(kind=c_int), value :: aero_rep_unique_id
197 !> Computational particle index
198 integer(kind=c_int), value :: particle_id
199 !> New number (m)
200 real(kind=c_double), value :: number_conc
202
203 !> Free an update data object
204 pure subroutine aero_rep_free_update_data(update_data) bind (c)
205 use iso_c_binding
206 !> Update data
207 type(c_ptr), value, intent(in) :: update_data
208 end subroutine aero_rep_free_update_data
209
210 end interface
211
212contains
213
214!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
215
216 !> Constructor for aero_rep_single_particle_t
217 function constructor() result (new_obj)
218
219 !> New aerosol representation
220 type(aero_rep_single_particle_t), pointer :: new_obj
221
222 allocate(new_obj)
223
224 end function constructor
225
226!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
227
228 !> Initialize the aerosol representation data, validating component data and
229 !! loading any required information from the \c
230 !! aero_rep_data_t::property_set. This routine should be called once for
231 !! each aerosol representation at the beginning of a model run after all
232 !! the input files have been read in. It ensures all data required during
233 !! the model run are included in the condensed data arrays.
234 subroutine initialize(this, aero_phase_set, spec_state_id)
235
236 !> Aerosol representation data
237 class(aero_rep_single_particle_t), intent(inout) :: this
238 !> The set of aerosol phases
239 type(aero_phase_data_ptr), pointer, intent(in) :: aero_phase_set(:)
240 !> Beginning state id for this aerosol representation in the model species
241 !! state array
242 integer(kind=i_kind), intent(in) :: spec_state_id
243
244 ! Unordered layer names (only used during initialization)
245 type(string_t), allocatable :: layer_names_unordered(:)
246 ! Cover names (only used during initialization)
247 type(string_t), allocatable :: cover_names_unordered(:)
248 ! Index in layer_names_unordered for each layer from inner- to outer-most
249 ! layer
250 integer(kind=i_kind), allocatable :: ordered_layer_id(:)
251
252 character(len=:), allocatable :: key_name, layer_name, layer_covers, &
253 phase_name
254 type(property_t), pointer :: layers, layer
255 type(property_ptr), allocatable :: phases(:)
256 integer(kind=i_kind) :: i_particle, i_phase, i_layer, i_aero, curr_id
257 integer(kind=i_kind) :: i_cover, j_phase, j_layer, i_map, curr_phase
258 integer(kind=i_kind) :: num_phases, num_int_param, num_float_param, &
259 num_particles
260 logical :: found
261
262 ! Get the maximum number of computational particles
263 key_name = "maximum computational particles"
264 call assert_msg(331697425, &
265 this%property_set%get_int(key_name, num_particles), &
266 "Missing maximum number of computational particles")
267
268 ! Get the set of layers
269 key_name = "layers"
270 call assert_msg(314292954, &
271 this%property_set%get_property_t(key_name, layers), &
272 "Missing layers for single-particle aerosol "// &
273 "representation '"//this%rep_name//"'")
274 call assert_msg(168669831, layers%size() .gt. 0, &
275 "No Layers specified for single-particle layer "// &
276 "aerosol representation '"//this%rep_name//"'")
277
278 ! Allocate space for the working arrays
279 allocate(phases(layers%size()))
280 allocate(cover_names_unordered(layers%size()))
281 allocate(layer_names_unordered(layers%size()))
282
283 ! Loop through the layers, adding names and counting the spaces needed
284 ! on the condensed data arrays, and counting the total phases instances
285 num_phases = 0
286 call layers%iter_reset()
287 do i_layer = 1, layers%size()
288
289 ! Get the layer properties
290 call assert_msg(303808978, layers%get_property_t(val=layer), &
291 "Invalid structure for layer '"// &
292 layer_names_unordered(i_layer)%string// &
293 "' in single-particle layer representation '"// &
294 this%rep_name//"'")
295
296 ! Get the layer name
297 key_name = "name"
298 call assert_msg(364496472, layer%get_string(key_name, layer_name), &
299 "Missing layer name in single-particle layer aerosol "// &
300 "representation '"//this%rep_name//"'")
301 layer_names_unordered(i_layer)%string = layer_name
302
303 ! Get the cover name
304 key_name = "covers"
305 call assert_msg(350939595, layer%get_string(key_name, layer_covers), &
306 "Missing cover name in layer'"// &
307 layer_names_unordered(i_layer)%string// &
308 "' in single-particle layer aerosol representation '"// &
309 this%rep_name//"'")
310 cover_names_unordered(i_layer)%string = layer_covers
311
312 ! Get the set of phases
313 key_name = "phases"
314 call assert_msg(647756433, &
315 layer%get_property_t(key_name, phases(i_layer)%val_), &
316 "Missing phases for layer '"// &
317 layer_names_unordered(i_layer)%string// &
318 "' in single-particle layer aerosol representation '"// &
319 this%rep_name//"'")
320
321 ! Add the phases to the counter
322 call assert_msg(002679882, phases(i_layer)%val_%size().gt.0, &
323 "No phases specified for layer '"// &
324 layer_names_unordered(i_layer)%string// &
325 "' in single-particle layer aerosol representation '"// &
326 this%rep_name//"'")
327
328 ! add to running total of phase count
329 num_phases = num_phases + phases(i_layer)%val_%size()
330
331 call layers%iter_next()
332 end do
333
334 ! get the map of layer names after reordering
335 ordered_layer_id = ordered_layer_ids(layer_names_unordered, &
336 cover_names_unordered)
337
338 ! set the layer names
339 allocate(this%layer_names_(size(ordered_layer_id)))
340 this%layer_names_(:) = layer_names_unordered(ordered_layer_id(:))
341
342 ! Allocate condensed data arrays
343 num_int_param = num_int_prop_ + 2 * layers%size() + 3 * num_phases
344 num_float_param = num_real_prop_
345 allocate(this%condensed_data_int(num_int_param))
346 allocate(this%condensed_data_real(num_float_param))
347 this%condensed_data_int(:) = int(0, kind=i_kind)
348 this%condensed_data_real(:) = real(0.0, kind=dp)
349
350 ! Save space for the environment-dependent parameters
351 this%num_env_params = num_env_param_per_particle_ * num_particles
352
353 ! Save representation dimensions
354 num_layers_ = layers%size()
355 max_particles_ = num_particles
356
357 ! validate phase names, assign aero_phase pointers for each phase in
358 ! each layer in each particle, and set PHASE_STATE_ID and
359 ! PHASE_MODEL_DATA_ID for each phase
360 allocate(this%aero_phase(num_phases * num_particles))
361 allocate(this%aero_phase_is_at_surface(num_phases * num_particles))
362 curr_phase = 1
363 do i_layer = 1, size(ordered_layer_id)
364 j_layer = ordered_layer_id(i_layer)
365
366 ! Set the starting and ending indices for the phases in this layer
367 layer_phase_start_(i_layer) = curr_phase
368 layer_phase_end_(i_layer) = curr_phase + phases(j_layer)%val_%size() - 1
369
370 curr_phase = curr_phase + phases(j_layer)%val_%size()
371 end do
372
373 curr_id = spec_state_id
374 this%state_id_start = spec_state_id
375 curr_phase = 1
376 do i_layer = 1, size(ordered_layer_id)
377 j_layer = ordered_layer_id(i_layer)
378
379 ! Loop through the phases and make sure they exist
380 call phases(j_layer)%val_%iter_reset()
381 do i_phase = 1, phases(j_layer)%val_%size()
382
383 ! Get the phase name
384 call assert_msg(566480284, &
385 phases(j_layer)%val_%get_string(val=phase_name), &
386 "Non-string phase name for layer '"// &
387 layer_names_unordered(j_layer)%string// &
388 "' in single-particle layer aerosol representation '"// &
389 this%rep_name//"'")
390
391 ! find phase and set pointer and indices
392 found = .false.
393 do j_phase = 1, size(aero_phase_set)
394 if (aero_phase_set(j_phase)%val%name() .eq. phase_name) then
395 found = .true.
396 do i_particle = 0, num_particles-1
397 this%aero_phase(i_particle*num_phases + curr_phase) = &
398 aero_phase_set(j_phase)
399 if (i_layer .eq. num_layers_) then
400 this%aero_phase_is_at_surface(i_particle*num_phases + curr_phase) = &
401 .true.
402 else
403 this%aero_phase_is_at_surface(i_particle*num_phases + curr_phase) = &
404 .false.
405 end if
406 end do
407 phase_state_id_(i_layer,i_phase) = curr_id
408 phase_model_data_id_(i_layer,i_phase) = j_phase
409 curr_id = curr_id + aero_phase_set(j_phase)%val%size()
410 curr_phase = curr_phase + 1
411 exit
412 end if
413 end do
414
415 call phases(j_layer)%val_%iter_next()
416 end do
417 end do
418 particle_state_size_ = curr_id - spec_state_id
419
420 ! Initialize the aerosol representation id
421 aero_rep_id_ = -1
422
423 ! Set the unique names for the chemical species
424 this%unique_names_ = this%unique_names( )
425
426 end subroutine initialize
427
428!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
429
430 !> Returns the maximum nunmber of computational particles
431 integer(kind=i_kind) function maximum_computational_particles(this)
432
433 !> Aerosol representation data
434 class(aero_rep_single_particle_t), intent(in) :: this
435
436 maximum_computational_particles = max_particles_
437
439
440!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
441
442 !> Get the size of the section of the
443 !! \c camp_camp_state::camp_state_t::state_var array required for this
444 !! aerosol representation.
445 !!
446 !! For a single particle representation, the size will correspond to the
447 !! the sum of the sizes of a single instance of each aerosol phase
448 !! provided to \c aero_rep_single_particle::initialize()
449 function get_size(this) result (state_size)
450
451 !> Size on the state array
452 integer(kind=i_kind) :: state_size
453 !> Aerosol representation data
454 class(aero_rep_single_particle_t), intent(in) :: this
455
456 state_size = max_particles_ * particle_state_size_
457
458 end function get_size
459
460!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
461
462 !> Get the number of state variables per-particle
463 !!
464 !! Calling functions can assume each particle has the same size on the
465 !! state array, and that individual particle states are contiguous and
466 !! arranged sequentially
467 function per_particle_size(this) result(state_size)
468
469 !> Size on the state array per particle
470 integer(kind=i_kind) :: state_size
471 !> Aerosol representation data
472 class(aero_rep_single_particle_t), intent(in) :: this
473
474 state_size = particle_state_size_
475
476 end function per_particle_size
477
478!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
479
480 !> Get a list of unique names for each element on the
481 !! \c camp_camp_state::camp_state_t::state_var array for this aerosol
482 !! representation. The list may be restricted to a particular phase and/or
483 !! aerosol species by including the phase_name and spec_name arguments.
484 !!
485 !! For a single particle representation, the unique names will be a 'P'
486 !! followed by the computational particle number, a '.', the phase name,
487 !! another '.', and the species name.
488 function unique_names(this, phase_name, tracer_type, spec_name)
489
490 use camp_util, only : integer_to_string
491 !> List of unique names
492 type(string_t), allocatable :: unique_names(:)
493 !> Aerosol representation data
494 class(aero_rep_single_particle_t), intent(in) :: this
495 !> Aerosol phase name
496 character(len=*), optional, intent(in) :: phase_name
497 !> Aerosol-phase species tracer type
498 integer(kind=i_kind), optional, intent(in) :: tracer_type
499 !> Aerosol-phase species name
500 character(len=*), optional, intent(in) :: spec_name
501
502 integer :: i_particle, i_layer, i_phase, i_spec, j_spec
503 integer :: num_spec, curr_tracer_type
504 type(string_t), allocatable :: spec_names(:)
505
506 ! copy saved unique names when available and no filters are included
507 if (.not. present(phase_name) .and. &
508 .not. present(tracer_type) .and. &
509 .not. present(spec_name) .and. &
510 allocated(this%unique_names_)) then
511 unique_names = this%unique_names_
512 return
513 end if
514
515 ! count the number of unique names
516 num_spec = 0
517 do i_phase = 1, size(this%aero_phase)
518 if (present(phase_name)) then
519 if(phase_name .ne. this%aero_phase(i_phase)%val%name()) cycle
520 end if
521 if (present(spec_name) .or. present(tracer_type)) then
522 spec_names = this%aero_phase(i_phase)%val%get_species_names()
523 do j_spec = 1, size(spec_names)
524 if (present(spec_name)) then
525 if (spec_name .ne. spec_names(j_spec)%string) cycle
526 end if
527 if (present(tracer_type)) then
528 curr_tracer_type = &
529 this%aero_phase(i_phase)%val%get_species_type( &
530 spec_names(j_spec)%string)
531 if (tracer_type .ne. curr_tracer_type) cycle
532 end if
533 num_spec = num_spec + 1
534 end do
535 deallocate(spec_names)
536 else
537 num_spec = num_spec + this%aero_phase(i_phase)%val%size()
538 end if
539 end do
540
541 ! allocate space for the unique names and assign them
542 num_spec = num_spec / max_particles_ ! we need per-particle value for indexing
543 allocate(unique_names(num_spec*max_particles_))
544 i_spec = 1
545 do i_layer = 1, num_layers_
546 do i_phase = layer_phase_start_(i_layer), layer_phase_end_(i_layer)
547 if (present(phase_name)) then
548 if(phase_name .ne. this%aero_phase(i_phase)%val%name()) cycle
549 end if
550 spec_names = this%aero_phase(i_phase)%val%get_species_names()
551 do j_spec = 1, this%aero_phase(i_phase)%val%size()
552 if (present(spec_name)) then
553 if (spec_name .ne. spec_names(j_spec)%string) cycle
554 end if
555 if (present(tracer_type)) then
556 curr_tracer_type = &
557 this%aero_phase(i_phase)%val%get_species_type( &
558 spec_names(j_spec)%string)
559 if (tracer_type .ne. curr_tracer_type) cycle
560 end if
561 do i_particle = 1, max_particles_
562 unique_names((i_particle-1)*num_spec+i_spec)%string = 'P'// &
563 trim(integer_to_string(i_particle))//"."// &
564 this%layer_names_(i_layer)%string//"."// &
565 this%aero_phase(i_phase)%val%name()//"."// &
566 spec_names(j_spec)%string
567 end do
568 i_spec = i_spec + 1
569 end do
570 deallocate(spec_names)
571 end do
572 end do
573
574 end function unique_names
575
576!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
577
578 !> Get a species id on the \c camp_camp_state::camp_state_t::state_var
579 !! array by its unique name. These are unique ids for each element on the
580 !! state array for this \ref camp_aero_rep "aerosol representation" and
581 !! are numbered:
582 !!
583 !! \f[x_u \in x_f ... (x_f+n-1)\f]
584 !!
585 !! where \f$x_u\f$ is the id of the element corresponding to the species
586 !! with unique name \f$u\f$ on the \c
587 !! camp_camp_state::camp_state_t::state_var array, \f$x_f\f$ is the index
588 !! of the first element for this aerosol representation on the state array
589 !! and \f$n\f$ is the total number of variables on the state array from
590 !! this aerosol representation.
591 function spec_state_id(this, unique_name) result (spec_id)
592
593 !> Species state id
594 integer(kind=i_kind) :: spec_id
595 !> Aerosol representation data
596 class(aero_rep_single_particle_t), intent(in) :: this
597 !> Unique name
598 character(len=*), intent(in) :: unique_name
599
600 integer(kind=i_kind) :: i_spec
601
602 spec_id = this%state_id_start
603 do i_spec = 1, size(this%unique_names_)
604 if (this%unique_names_(i_spec)%string .eq. unique_name) then
605 return
606 end if
607 spec_id = spec_id + 1
608 end do
609 call die_msg( 449087541, "Cannot find species '"//unique_name//"'" )
610
611 end function spec_state_id
612
613!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
614
615 !> Get the non-unique name of a species in this aerosol representation by
616 !! id.
617 function spec_name(this, unique_name)
618
619 !> Chemical species name
620 character(len=:), allocatable :: spec_name
621 !> Aerosol representation data
622 class(aero_rep_single_particle_t), intent(in) :: this
623 !> Unique name of the species in this aerosol representation
624 character(len=*), intent(in) :: unique_name
625
626 type(string_t) :: l_unique_name
627 type(string_t), allocatable :: substrs(:)
628
629 l_unique_name%string = unique_name
630 substrs = l_unique_name%split(".")
631 call assert(407537518, size( substrs ) .eq. 4 )
632 spec_name = substrs(4)%string
633
634 end function spec_name
635
636!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
637
638 !> Get the number of instances of a specified aerosol phase. In the single
639 !! particle representation with layers, a phase can exist in multiple layers
640 !! in one particle.
641 integer(kind=i_kind) function num_phase_instances(this, phase_name)
642
643 !> Aerosol representation data
644 class(aero_rep_single_particle_t), intent(in) :: this
645 !> Aerosol phase name
646 character(len=*), intent(in) :: phase_name
647
648 integer(kind=i_kind) :: i_phase, i_layer, phase_index
649
651 phase_index = 0
652 do i_layer = 1, num_layers_
653 do i_phase = layer_phase_start_(i_layer), layer_phase_end_(i_layer)
654 if (this%aero_phase(i_phase)%val%name() .eq. phase_name) then
655 phase_index = phase_index + 1
656 end if
657 end do
658 end do
659 num_phase_instances = phase_index * max_particles_
660
661 end function num_phase_instances
662
663!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
664
665 !> Get the number of Jacobian elements used in calculations of aerosol mass,
666 !! volume, number, etc. for a particular phase
667 function num_jac_elem(this, phase_id)
668
669 !> Number of Jacobian elements used
670 integer(kind=i_kind) :: num_jac_elem
671 !> Aerosol respresentation data
672 class(aero_rep_single_particle_t), intent(in) :: this
673 !> Aerosol phase id
674 integer(kind=i_kind), intent(in) :: phase_id
675
676 integer(kind=i_kind) :: i_phase
677
678 call assert_msg(927040495, phase_id .ge. 1 .and. &
679 phase_id .le. size( this%aero_phase ), &
680 "Aerosol phase index out of range. Got "// &
681 trim( integer_to_string( phase_id ) )//", expected 1:"// &
682 trim( integer_to_string( size( this%aero_phase ) ) ) )
683 num_jac_elem = 0
684 do i_phase = 1, total_num_phases_
686 this%aero_phase(i_phase)%val%num_jac_elem()
687 end do
688
689 end function num_jac_elem
690
691!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
692
693 !> Returns the number of layers
694 integer function num_layers(this)
695
696 !> Aerosol representation data
697 class(aero_rep_single_particle_t), intent(in) :: this
698
699 num_layers = num_layers_
700
701 end function num_layers
702
703!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
704
705 !> Returns the number of phases in a layer or overall
706 integer function num_phases(this, layer)
707
708 !> Aerosol representation data
709 class(aero_rep_single_particle_t), intent(in) :: this
710 !> Layer id
711 integer, optional, intent(in) :: layer
712
713 if (present(layer)) then
714 num_phases = layer_phase_end_(layer) - layer_phase_start_(layer) + 1
715 else
716 num_phases = layer_phase_end_(num_layers_) - layer_phase_start_(1) + 1
717 end if
718
719 end function num_phases
720
721!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
722
723 !> Returns the number of state variables for a layer and phase
724 integer function phase_state_size(this, layer, phase)
725
726 use camp_util, only : die_msg
727
728 !> Aerosol representation data
729 class(aero_rep_single_particle_t), intent(in) :: this
730 !> Layer id
731 integer, optional, intent(in) :: layer
732 !> Phase id
733 integer, optional, intent(in) :: phase
734
735 if (present(layer) .and. present(phase)) then
736 if (layer .eq. num_layers_ .and. phase .eq. num_phases_(layer)) then
737 phase_state_size = phase_state_id_(1,1) + particle_state_size_ - &
738 phase_state_id_(layer, phase)
739 else if (phase .eq. num_phases_(layer)) then
740 phase_state_size = phase_state_id_(layer+1, 1) - &
741 phase_state_id_(layer, phase)
742 else
743 phase_state_size = phase_state_id_(layer, phase+1) - &
744 phase_state_id_(layer, phase)
745 end if
746 else if (present(layer)) then
747 if (layer .eq. num_layers_) then
748 phase_state_size = phase_state_id_(1,1) + particle_state_size_ - &
749 phase_state_id_(layer, 1)
750 else
751 phase_state_size = phase_state_id_(layer+1, 1) - &
752 phase_state_id_(layer, 1)
753 end if
754 else if (present(phase)) then
755 call die_msg(917793122, "Must specify layer if including phase is "// &
756 "state size request")
757 else
758 phase_state_size = particle_state_size_
759 end if
760
761 end function phase_state_size
762
763!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
764
765 !> Finalize the aerosol representation
766 elemental subroutine finalize(this)
767
768 !> Aerosol representation data
769 type(aero_rep_single_particle_t), intent(inout) :: this
770
771 if (allocated(this%rep_name)) deallocate(this%rep_name)
772 if (allocated(this%aero_phase)) then
773 ! The core will deallocate the aerosol phases
774 call this%aero_phase(:)%dereference()
775 deallocate(this%aero_phase)
776 end if
777 if (allocated(this%unique_names_)) deallocate(this%unique_names_)
778 if (allocated(this%layer_names_)) deallocate(this%layer_names_)
779 if (associated(this%property_set)) deallocate(this%property_set)
780 if (allocated(this%condensed_data_real)) &
781 deallocate(this%condensed_data_real)
782 if (allocated(this%condensed_data_int)) &
783 deallocate(this%condensed_data_int)
784
785 end subroutine finalize
786
787!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
788
789 !> Initialize an update data object
790 subroutine update_data_init_number(this, update_data, aero_rep_type)
791
792 use camp_rand, only : generate_int_id
793
794 !> Aerosol representation to update
795 class(aero_rep_single_particle_t), intent(inout) :: this
796 !> Update data object
797 class(aero_rep_update_data_single_particle_number_t), intent(out) :: &
798 update_data
799 !> Aerosol representaiton id
800 integer(kind=i_kind), intent(in) :: aero_rep_type
801
802 ! If an aerosol representation id has not been generated, do it now
803 if (aero_rep_id_.eq.-1) then
804 aero_rep_id_ = generate_int_id()
805 end if
806
807 update_data%aero_rep_unique_id = aero_rep_id_
808 update_data%maximum_computational_particles = &
809 this%maximum_computational_particles( )
810 update_data%aero_rep_type = int(aero_rep_type, kind=c_int)
811 update_data%update_data = &
813 update_data%is_malloced = .true.
814
815 end subroutine update_data_init_number
816
817!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
818
819 !> Order layer array from inner most layer to outermost
820 function ordered_layer_ids(layer_names_unordered, cover_names_unordered)
821
822 !> Layer names in original order
823 type(string_t), intent(in) :: layer_names_unordered(:)
824 !> Name of "covered" layer for each layer in layer_name_unordered
825 type(string_t), intent(in) :: cover_names_unordered(:)
826 !> Index of name in layer_name_unordered for each layer after reordering
827 integer, allocatable :: ordered_layer_ids(:)
828
829 integer(kind=i_kind) :: i_layer, j_layer, i_cover
830
831 ! Ensure layer names do not repeat
832 do i_layer = 1, size(layer_names_unordered)
833 do j_layer = 1, size(layer_names_unordered)
834 if (i_layer .eq. j_layer) cycle
835 call assert_msg(781626922, layer_names_unordered(i_layer)%string .ne. &
836 layer_names_unordered(j_layer)%string, &
837 "Duplicate layer name in single particle "// &
838 "representation: '"// &
839 trim(layer_names_unordered(i_layer)%string)//"'")
840 end do
841 end do
842
843 allocate(ordered_layer_ids(size(layer_names_unordered)))
844
845 ! Search for innermost layer with cover set to 'none'
846 do i_layer = 1, size(layer_names_unordered)
847 if (cover_names_unordered(i_layer)%string == "none") then
848 ordered_layer_ids(1) = i_layer
849 end if
850 end do
851
852 ! Assign each layer working outwards from center of particle
853 do i_cover = 2, size(ordered_layer_ids)
854 do i_layer = 1, size(layer_names_unordered)
855 if (layer_names_unordered(ordered_layer_ids(i_cover-1))%string &
856 .eq. cover_names_unordered(i_layer)%string) then
857 ordered_layer_ids(i_cover) = i_layer
858 exit
859 end if
860 end do
861 end do
862
863 end function ordered_layer_ids
864
865!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
866
867 !> Set packed update data for particle number (#/m3) for a particular
868 !! computational particle.
869 subroutine update_data_set_number__n_m3(this, particle_id, number_conc)
870
871 !> Update data
872 class(aero_rep_update_data_single_particle_number_t), intent(inout) :: &
873 this
874 !> Computational particle index
875 integer(kind=i_kind), intent(in) :: particle_id
876 !> Updated number
877 real(kind=dp), intent(in) :: number_conc
878
879 call assert_msg(611967802, this%is_malloced, &
880 "Trying to set number of uninitialized update object.")
881 call assert_msg(689085496, particle_id .ge. 1 .and. &
882 particle_id .le. this%maximum_computational_particles, &
883 "Invalid computational particle index: "// &
884 trim(integer_to_string(particle_id)))
886 this%get_data(), this%aero_rep_unique_id, particle_id-1, &
887 number_conc)
888
889 end subroutine update_data_set_number__n_m3
890
891!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
892
893 !> Determine the size of a binary required to pack the reaction data
894 integer(kind=i_kind) function internal_pack_size_number(this, comm) &
895 result(pack_size)
896
897 !> Aerosol representation update data
898 class(aero_rep_update_data_single_particle_number_t), intent(in) :: this
899 !> MPI communicator
900 integer, intent(in) :: comm
901
902 pack_size = &
903 camp_mpi_pack_size_logical(this%is_malloced, comm) + &
904 camp_mpi_pack_size_integer(this%maximum_computational_particles, comm) + &
905 camp_mpi_pack_size_integer(this%aero_rep_unique_id, comm)
906
907 end function internal_pack_size_number
908
909!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
910
911 !> Pack the given value to the buffer, advancing position
912 subroutine internal_bin_pack_number(this, buffer, pos, comm)
913
914 !> Aerosol representation update data
915 class(aero_rep_update_data_single_particle_number_t), intent(in) :: this
916 !> Memory buffer
917 character, intent(inout) :: buffer(:)
918 !> Current buffer position
919 integer, intent(inout) :: pos
920 !> MPI communicator
921 integer, intent(in) :: comm
922
923#ifdef CAMP_USE_MPI
924 integer :: prev_position
925
926 prev_position = pos
927 call camp_mpi_pack_logical(buffer, pos, this%is_malloced, comm)
928 call camp_mpi_pack_integer(buffer, pos, &
929 this%maximum_computational_particles, comm)
930 call camp_mpi_pack_integer(buffer, pos, this%aero_rep_unique_id, comm)
931 call assert(411585487, &
932 pos - prev_position <= this%pack_size(comm))
933#endif
934
935 end subroutine internal_bin_pack_number
936
937!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
938
939 !> Unpack the given value from the buffer, advancing position
940 subroutine internal_bin_unpack_number(this, buffer, pos, comm)
941
942 !> Aerosol representation update data
943 class(aero_rep_update_data_single_particle_number_t), intent(inout) :: this
944 !> Memory buffer
945 character, intent(inout) :: buffer(:)
946 !> Current buffer position
947 integer, intent(inout) :: pos
948 !> MPI communicator
949 integer, intent(in) :: comm
950
951#ifdef CAMP_USE_MPI
952 integer :: prev_position
953
954 prev_position = pos
955 call camp_mpi_unpack_logical(buffer, pos, this%is_malloced, comm)
956 call camp_mpi_unpack_integer(buffer, pos, &
957 this%maximum_computational_particles, comm)
958 call camp_mpi_unpack_integer(buffer, pos, this%aero_rep_unique_id, comm)
959 call assert(351557153, &
960 pos - prev_position <= this%pack_size(comm))
962#endif
963
964 end subroutine internal_bin_unpack_number
965
966!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
967
968 !> Finalize a number update data object
969 elemental subroutine update_data_number_finalize(this)
970
971 !> Update data object to free
972 type(aero_rep_update_data_single_particle_number_t), intent(inout) :: this
973
974 if (this%is_malloced) call aero_rep_free_update_data(this%update_data)
975
976 end subroutine update_data_number_finalize
977
978!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
979
Get the size of the section of the camp_camp_state::camp_state_t::state_var array required for this a...
Initialize the aerosol representation data, validating component data and loading any required inform...
Extending-type binary pack function (Internal use only)
Extending-type binary unpack function (Internal use only)
Extending-type binary pack size (internal use only)
Get the number of Jacobian elements used in calculations of aerosol mass, volume, number,...
Get the number of instances of a specified aerosol phase.
Get the non-unique name of a chemical species by its unique name.
Get a species id on the camp_camp_state::camp_state_t::state_var array by unique name....
Get a list of unique names for each element on the camp_camp_state::camp_state_t::state_var array for...
Interface for to_string functions.
Definition util.F90:32
The abstract aero_phase_data_t structure and associated subroutines.
elemental subroutine finalize(this)
Finalize the aerosol phase data.
type(aero_phase_data_t) function, pointer constructor(phase_name, init_size)
Constructor for aero_phase_data_t.
integer(kind=i_kind) function pack_size(this, comm)
Determine the size of a binary required to pack the aerosol representation data.
The abstract aero_rep_data_t structure and associated subroutines.
The 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.
elemental subroutine update_data_number_finalize(this)
Finalize a number update data object.
integer function num_layers(this)
Returns the number of layers.
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_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.
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:38