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
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 ! Loop through the phases and make sure they exist
379 call phases(j_layer)%val_%iter_reset()
380 do i_phase = 1, phases(j_layer)%val_%size()
381
382 ! Get the phase name
383 call assert_msg(566480284, &
384 phases(j_layer)%val_%get_string(val=phase_name), &
385 "Non-string phase name for layer '"// &
386 layer_names_unordered(j_layer)%string// &
387 "' in single-particle layer aerosol representation '"// &
388 this%rep_name//"'")
389 ! find phase and set pointer and indices
390 found = .false.
391 do j_phase = 1, size(aero_phase_set)
392 if (aero_phase_set(j_phase)%val%name() .eq. phase_name) then
393 found = .true.
394 do i_particle = 0, num_particles-1
395 this%aero_phase(i_particle*num_phases + curr_phase) = &
396 aero_phase_set(j_phase)
397 if (i_layer .eq. num_layers_) then
398 this%aero_phase_is_at_surface(i_particle*num_phases + curr_phase) = &
399 .true.
400 else
401 this%aero_phase_is_at_surface(i_particle*num_phases + curr_phase) = &
402 .false.
403 end if
404 end do
405 phase_state_id_(i_layer,i_phase) = curr_id
406 phase_model_data_id_(i_layer,i_phase) = j_phase
407 curr_id = curr_id + aero_phase_set(j_phase)%val%size()
408 curr_phase = curr_phase + 1
409 exit
410 end if
411 end do
412 call assert(373124707, found)
413
414 call phases(j_layer)%val_%iter_next()
415 end do
416 end do
417 particle_state_size_ = curr_id - spec_state_id
418
419 ! Initialize the aerosol representation id
420 aero_rep_id_ = -1
421
422 ! Set the unique names for the chemical species
423 this%unique_names_ = this%unique_names( )
424
425 end subroutine initialize
426
427!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
428
429 !> Returns the maximum nunmber of computational particles
430 integer(kind=i_kind) function maximum_computational_particles(this)
431
432 !> Aerosol representation data
433 class(aero_rep_single_particle_t), intent(in) :: this
434
435 maximum_computational_particles = max_particles_
436
438
439!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
440
441 !> Get the size of the section of the
442 !! \c camp_camp_state::camp_state_t::state_var array required for this
443 !! aerosol representation.
444 !!
445 !! For a single particle representation, the size will correspond to the
446 !! the sum of the sizes of a single instance of each aerosol phase
447 !! provided to \c aero_rep_single_particle::initialize()
448 function get_size(this) result (state_size)
449
450 !> Size on the state array
451 integer(kind=i_kind) :: state_size
452 !> Aerosol representation data
453 class(aero_rep_single_particle_t), intent(in) :: this
454
455 state_size = max_particles_ * particle_state_size_
456
457 end function get_size
458
459!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
460
461 !> Get the number of state variables per-particle
462 !!
463 !! Calling functions can assume each particle has the same size on the
464 !! state array, and that individual particle states are contiguous and
465 !! arranged sequentially
466 function per_particle_size(this) result(state_size)
467
468 !> Size on the state array per particle
469 integer(kind=i_kind) :: state_size
470 !> Aerosol representation data
471 class(aero_rep_single_particle_t), intent(in) :: this
472
473 state_size = particle_state_size_
474
475 end function per_particle_size
476
477!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
478
479 !> Get a list of unique names for each element on the
480 !! \c camp_camp_state::camp_state_t::state_var array for this aerosol
481 !! representation. The list may be restricted to a particular phase and/or
482 !! aerosol species by including the phase_name and spec_name arguments.
483 !!
484 !! For a single particle representation, the unique names will be a 'P'
485 !! followed by the computational particle number, a '.', the phase name,
486 !! another '.', and the species name.
487 function unique_names(this, phase_name, tracer_type, spec_name, &
488 phase_is_at_surface)
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 !> Indicates if aerosol phase is at the surface of particle
502 logical, optional, intent(in) :: phase_is_at_surface
503
504 integer :: i_particle, i_layer, i_phase, i_spec, j_spec
505 integer :: num_spec, curr_tracer_type
506 type(string_t), allocatable :: spec_names(:)
507
508 ! copy saved unique names when available and no filters are included
509 if (.not. present(phase_name) .and. &
510 .not. present(tracer_type) .and. &
511 .not. present(spec_name) .and. &
512 allocated(this%unique_names_)) then
513 unique_names = this%unique_names_
514 return
515 end if
516
517 ! count the number of unique names
518 num_spec = 0
519 do i_phase = 1, size(this%aero_phase)
520 if (present(phase_name)) then
521 if(phase_name .ne. this%aero_phase(i_phase)%val%name()) cycle
522 end if
523 if (present(phase_is_at_surface)) then
524 if (phase_is_at_surface .neqv. &
525 this%aero_phase_is_at_surface(i_phase)) cycle
526 end if
527 if (present(spec_name) .or. present(tracer_type)) then
528 spec_names = this%aero_phase(i_phase)%val%get_species_names()
529 do j_spec = 1, size(spec_names)
530 if (present(spec_name)) then
531 if (spec_name .ne. spec_names(j_spec)%string) cycle
532 end if
533 if (present(tracer_type)) then
534 curr_tracer_type = &
535 this%aero_phase(i_phase)%val%get_species_type( &
536 spec_names(j_spec)%string)
537 if (tracer_type .ne. curr_tracer_type) cycle
538 end if
539 num_spec = num_spec + 1
540 end do
541 deallocate(spec_names)
542 else
543 num_spec = num_spec + this%aero_phase(i_phase)%val%size()
544 end if
545 end do
546
547 ! allocate space for the unique names and assign them
548 num_spec = num_spec / max_particles_ ! we need per-particle value for indexing
549 allocate(unique_names(num_spec*max_particles_))
550 i_spec = 1
551 do i_layer = 1, num_layers_
552 do i_phase = layer_phase_start_(i_layer), layer_phase_end_(i_layer)
553 if (present(phase_name)) then
554 if(phase_name .ne. this%aero_phase(i_phase)%val%name()) cycle
555 end if
556 if (present(phase_is_at_surface)) then
557 if (phase_is_at_surface .neqv. &
558 this%aero_phase_is_at_surface(i_phase)) cycle
559 end if
560 spec_names = this%aero_phase(i_phase)%val%get_species_names()
561 do j_spec = 1, this%aero_phase(i_phase)%val%size()
562 if (present(spec_name)) then
563 if (spec_name .ne. spec_names(j_spec)%string) cycle
564 end if
565 if (present(tracer_type)) then
566 curr_tracer_type = &
567 this%aero_phase(i_phase)%val%get_species_type( &
568 spec_names(j_spec)%string)
569 if (tracer_type .ne. curr_tracer_type) cycle
570 end if
571 do i_particle = 1, max_particles_
572 unique_names((i_particle-1)*num_spec+i_spec)%string = 'P'// &
573 trim(integer_to_string(i_particle))//"."// &
574 this%layer_names_(i_layer)%string//"."// &
575 this%aero_phase(i_phase)%val%name()//"."// &
576 spec_names(j_spec)%string
577 end do
578 i_spec = i_spec + 1
579 end do
580 deallocate(spec_names)
581 end do
582 end do
583
584 end function unique_names
585
586!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
587
588 !> Get a species id on the \c camp_camp_state::camp_state_t::state_var
589 !! array by its unique name. These are unique ids for each element on the
590 !! state array for this \ref camp_aero_rep "aerosol representation" and
591 !! are numbered:
592 !!
593 !! \f[x_u \in x_f ... (x_f+n-1)\f]
594 !!
595 !! where \f$x_u\f$ is the id of the element corresponding to the species
596 !! with unique name \f$u\f$ on the \c
597 !! camp_camp_state::camp_state_t::state_var array, \f$x_f\f$ is the index
598 !! of the first element for this aerosol representation on the state array
599 !! and \f$n\f$ is the total number of variables on the state array from
600 !! this aerosol representation.
601 function spec_state_id(this, unique_name) result (spec_id)
602
603 !> Species state id
604 integer(kind=i_kind) :: spec_id
605 !> Aerosol representation data
606 class(aero_rep_single_particle_t), intent(in) :: this
607 !> Unique name
608 character(len=*), intent(in) :: unique_name
609
610 integer(kind=i_kind) :: i_spec
611
612 spec_id = this%state_id_start
613 do i_spec = 1, size(this%unique_names_)
614 if (this%unique_names_(i_spec)%string .eq. unique_name) then
615 return
616 end if
617 spec_id = spec_id + 1
618 end do
619 call die_msg( 449087541, "Cannot find species '"//unique_name//"'" )
620
621 end function spec_state_id
622
623!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
624
625 !> Get the non-unique name of a species in this aerosol representation by
626 !! id.
627 function spec_name(this, unique_name)
628
629 !> Chemical species name
630 character(len=:), allocatable :: spec_name
631 !> Aerosol representation data
632 class(aero_rep_single_particle_t), intent(in) :: this
633 !> Unique name of the species in this aerosol representation
634 character(len=*), intent(in) :: unique_name
635
636 type(string_t) :: l_unique_name
637 type(string_t), allocatable :: substrs(:)
638
639 l_unique_name%string = unique_name
640 substrs = l_unique_name%split(".")
641 call assert(407537518, size( substrs ) .eq. 4 )
642 spec_name = substrs(4)%string
643
644 end function spec_name
645
646!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
647
648 !> Get the number of instances of a specified aerosol phase. In the single
649 !! particle representation with layers, a phase can exist in multiple layers
650 !! in one particle.
651 integer(kind=i_kind) function num_phase_instances(this, phase_name, &
652 is_at_surface)
653
654 !> Aerosol representation data
655 class(aero_rep_single_particle_t), intent(in) :: this
656 !> Aerosol phase name
657 character(len=*), intent(in) :: phase_name
658 !> Indicates if aerosol phase is at the surface of particle
659 logical, intent(in), optional :: is_at_surface
660
661 integer(kind=i_kind) :: i_phase, i_layer, phase_index
662
664 if (present(is_at_surface)) then
665 do i_layer = 1, num_layers_
666 do i_phase = layer_phase_start_(i_layer), layer_phase_end_(i_layer)
667 if (this%aero_phase(i_phase)%val%name() .eq. phase_name) then
668 if (this%aero_phase_is_at_surface(i_phase) .eqv. is_at_surface) then
670 end if
671 end if
672 end do
673 end do
674 else
675 do i_layer = 1, num_layers_
676 do i_phase = layer_phase_start_(i_layer), layer_phase_end_(i_layer)
677 if (this%aero_phase(i_phase)%val%name() .eq. phase_name) then
679 end if
680 end do
681 end do
682 end if
684
685 end function num_phase_instances
686
687!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
688
689 !> Get the number of Jacobian elements used in calculations of aerosol mass,
690 !! volume, number, etc. for a particular phase
691 function num_jac_elem(this, phase_id)
692
693 !> Number of Jacobian elements used
694 integer(kind=i_kind) :: num_jac_elem
695 !> Aerosol respresentation data
696 class(aero_rep_single_particle_t), intent(in) :: this
697 !> Aerosol phase id
698 integer(kind=i_kind), intent(in) :: phase_id
699
700 integer(kind=i_kind) :: i_phase
701
702 call assert_msg(927040495, phase_id .ge. 1 .and. &
703 phase_id .le. size( this%aero_phase ), &
704 "Aerosol phase index out of range. Got "// &
705 trim( integer_to_string( phase_id ) )//", expected 1:"// &
706 trim( integer_to_string( size( this%aero_phase ) ) ) )
707 num_jac_elem = 0
708 do i_phase = 1, total_num_phases_
710 this%aero_phase(i_phase)%val%num_jac_elem()
711 end do
712
713 end function num_jac_elem
714
715!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
716
717 !> Returns the number of layers
718 integer function num_layers(this)
719
720 !> Aerosol representation data
721 class(aero_rep_single_particle_t), intent(in) :: this
722
723 num_layers = num_layers_
724
725 end function num_layers
726
727!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
728
729 !> Returns the number of phases in a layer or overall
730 integer function num_phases(this, layer)
731
732 !> Aerosol representation data
733 class(aero_rep_single_particle_t), intent(in) :: this
734 !> Layer id
735 integer, optional, intent(in) :: layer
736
737 if (present(layer)) then
738 num_phases = layer_phase_end_(layer) - layer_phase_start_(layer) + 1
739 else
740 num_phases = layer_phase_end_(num_layers_) - layer_phase_start_(1) + 1
741 end if
742
743 end function num_phases
744
745!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
746
747 !> Returns the number of state variables for a layer and phase
748 integer function phase_state_size(this, layer, phase)
749
750 use camp_util, only : die_msg
751
752 !> Aerosol representation data
753 class(aero_rep_single_particle_t), intent(in) :: this
754 !> Layer id
755 integer, optional, intent(in) :: layer
756 !> Phase id
757 integer, optional, intent(in) :: phase
758
759 if (present(layer) .and. present(phase)) then
760 if (layer .eq. num_layers_ .and. phase .eq. num_phases_(layer)) then
761 phase_state_size = phase_state_id_(1,1) + particle_state_size_ - &
762 phase_state_id_(layer, phase)
763 else if (phase .eq. num_phases_(layer)) then
764 phase_state_size = phase_state_id_(layer+1, 1) - &
765 phase_state_id_(layer, phase)
766 else
767 phase_state_size = phase_state_id_(layer, phase+1) - &
768 phase_state_id_(layer, phase)
769 end if
770 else if (present(layer)) then
771 if (layer .eq. num_layers_) then
772 phase_state_size = phase_state_id_(1,1) + particle_state_size_ - &
773 phase_state_id_(layer, 1)
774 else
775 phase_state_size = phase_state_id_(layer+1, 1) - &
776 phase_state_id_(layer, 1)
777 end if
778 else if (present(phase)) then
779 call die_msg(917793122, "Must specify layer if including phase is "// &
780 "state size request")
781 else
782 phase_state_size = particle_state_size_
783 end if
784
785 end function phase_state_size
786
787!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
788
789 !> Finalize the aerosol representation
790 subroutine finalize(this)
791
792 !> Aerosol representation data
793 type(aero_rep_single_particle_t), intent(inout) :: this
794
795 if (allocated(this%rep_name)) deallocate(this%rep_name)
796 if (allocated(this%aero_phase)) then
797 ! The core will deallocate the aerosol phases
798 call this%aero_phase(:)%dereference()
799 deallocate(this%aero_phase)
800 end if
801 if (allocated(this%unique_names_)) deallocate(this%unique_names_)
802 if (allocated(this%layer_names_)) deallocate(this%layer_names_)
803 if (associated(this%property_set)) deallocate(this%property_set)
804 if (allocated(this%condensed_data_real)) &
805 deallocate(this%condensed_data_real)
806 if (allocated(this%condensed_data_int)) &
807 deallocate(this%condensed_data_int)
808
809 end subroutine finalize
810
811!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
812
813 !> Finalize the aerosol representation array
814 subroutine finalize_array(this_array)
815
816 !> Aerosol representation data
817 type(aero_rep_single_particle_t), intent(inout) :: this_array(:)
818
819 integer(kind=i_kind) :: i_aero
820
821 do i_aero = 1, size(this_array)
822 call finalize(this_array(i_aero))
823 end do
824
825 end subroutine finalize_array
826
827!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
828
829 !> Initialize an update data object
830 subroutine update_data_init_number(this, update_data, aero_rep_type)
831
832 use camp_rand, only : generate_int_id
833
834 !> Aerosol representation to update
835 class(aero_rep_single_particle_t), intent(inout) :: this
836 !> Update data object
837 class(aero_rep_update_data_single_particle_number_t), intent(out) :: &
838 update_data
839 !> Aerosol representaiton id
840 integer(kind=i_kind), intent(in) :: aero_rep_type
841
842 ! If an aerosol representation id has not been generated, do it now
843 if (aero_rep_id_.eq.-1) then
844 aero_rep_id_ = generate_int_id()
845 end if
846
847 update_data%aero_rep_unique_id = aero_rep_id_
848 update_data%maximum_computational_particles = &
849 this%maximum_computational_particles( )
850 update_data%aero_rep_type = int(aero_rep_type, kind=c_int)
851 update_data%update_data = &
853 update_data%is_malloced = .true.
854
855 end subroutine update_data_init_number
856
857!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
858
859 !> Order layer array from inner most layer to outermost
860 function ordered_layer_ids(layer_names_unordered, cover_names_unordered)
861
862 !> Layer names in original order
863 type(string_t), intent(in) :: layer_names_unordered(:)
864 !> Name of "covered" layer for each layer in layer_name_unordered
865 type(string_t), intent(in) :: cover_names_unordered(:)
866 !> Index of name in layer_name_unordered for each layer after reordering
867 integer, allocatable :: ordered_layer_ids(:)
868
869 integer(kind=i_kind) :: i_layer, j_layer, i_cover
870
871 ! Ensure layer names do not repeat
872 do i_layer = 1, size(layer_names_unordered)
873 do j_layer = 1, size(layer_names_unordered)
874 if (i_layer .eq. j_layer) cycle
875 call assert_msg(781626922, layer_names_unordered(i_layer)%string .ne. &
876 layer_names_unordered(j_layer)%string, &
877 "Duplicate layer name in single particle "// &
878 "representation: '"// &
879 trim(layer_names_unordered(i_layer)%string)//"'")
880 end do
881 end do
882
883 allocate(ordered_layer_ids(size(layer_names_unordered)))
884
885 ! Search for innermost layer with cover set to 'none'
886 do i_layer = 1, size(layer_names_unordered)
887 if (cover_names_unordered(i_layer)%string == "none") then
888 ordered_layer_ids(1) = i_layer
889 end if
890 end do
891
892 ! Assign each layer working outwards from center of particle
893 do i_cover = 2, size(ordered_layer_ids)
894 do i_layer = 1, size(layer_names_unordered)
895 if (layer_names_unordered(ordered_layer_ids(i_cover-1))%string &
896 .eq. cover_names_unordered(i_layer)%string) then
897 ordered_layer_ids(i_cover) = i_layer
898 exit
899 end if
900 end do
901 end do
902
903 end function ordered_layer_ids
904
905!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
906
907 !> Set packed update data for particle number (#/m3) for a particular
908 !! computational particle.
909 subroutine update_data_set_number__n_m3(this, particle_id, number_conc)
910
911 !> Update data
912 class(aero_rep_update_data_single_particle_number_t), intent(inout) :: &
913 this
914 !> Computational particle index
915 integer(kind=i_kind), intent(in) :: particle_id
916 !> Updated number
917 real(kind=dp), intent(in) :: number_conc
918
919 call assert_msg(611967802, this%is_malloced, &
920 "Trying to set number of uninitialized update object.")
921 call assert_msg(689085496, particle_id .ge. 1 .and. &
922 particle_id .le. this%maximum_computational_particles, &
923 "Invalid computational particle index: "// &
924 trim(integer_to_string(particle_id)))
926 this%get_data(), this%aero_rep_unique_id, particle_id-1, &
927 number_conc)
928
929 end subroutine update_data_set_number__n_m3
930
931!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
932
933 !> Determine the size of a binary required to pack the reaction data
934 integer(kind=i_kind) function internal_pack_size_number(this, comm) &
935 result(pack_size)
936
937 !> Aerosol representation update data
938 class(aero_rep_update_data_single_particle_number_t), intent(in) :: this
939 !> MPI communicator
940 integer, intent(in) :: comm
941
942 pack_size = &
943 camp_mpi_pack_size_logical(this%is_malloced, comm) + &
944 camp_mpi_pack_size_integer(this%maximum_computational_particles, comm) + &
945 camp_mpi_pack_size_integer(this%aero_rep_unique_id, comm)
946
947 end function internal_pack_size_number
948
949!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
950
951 !> Pack the given value to the buffer, advancing position
952 subroutine internal_bin_pack_number(this, buffer, pos, comm)
953
954 !> Aerosol representation update data
955 class(aero_rep_update_data_single_particle_number_t), intent(in) :: this
956 !> Memory buffer
957 character, intent(inout) :: buffer(:)
958 !> Current buffer position
959 integer, intent(inout) :: pos
960 !> MPI communicator
961 integer, intent(in) :: comm
962
963#ifdef CAMP_USE_MPI
964 integer :: prev_position
965
966 prev_position = pos
967 call camp_mpi_pack_logical(buffer, pos, this%is_malloced, comm)
968 call camp_mpi_pack_integer(buffer, pos, &
969 this%maximum_computational_particles, comm)
970 call camp_mpi_pack_integer(buffer, pos, this%aero_rep_unique_id, comm)
971 call assert(411585487, &
972 pos - prev_position <= this%pack_size(comm))
973#endif
974
975 end subroutine internal_bin_pack_number
976
977!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
978
979 !> Unpack the given value from the buffer, advancing position
980 subroutine internal_bin_unpack_number(this, buffer, pos, comm)
981
982 !> Aerosol representation update data
983 class(aero_rep_update_data_single_particle_number_t), intent(inout) :: this
984 !> Memory buffer
985 character, intent(inout) :: buffer(:)
986 !> Current buffer position
987 integer, intent(inout) :: pos
988 !> MPI communicator
989 integer, intent(in) :: comm
990
991#ifdef CAMP_USE_MPI
992 integer :: prev_position
993
994 prev_position = pos
995 call camp_mpi_unpack_logical(buffer, pos, this%is_malloced, comm)
996 call camp_mpi_unpack_integer(buffer, pos, &
997 this%maximum_computational_particles, comm)
998 call camp_mpi_unpack_integer(buffer, pos, this%aero_rep_unique_id, comm)
999 call assert(351557153, &
1000 pos - prev_position <= this%pack_size(comm))
1002#endif
1003
1004 end subroutine internal_bin_unpack_number
1005
1006!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1007
1008 !> Finalize a number update data object
1010
1011 !> Update data object to free
1012 type(aero_rep_update_data_single_particle_number_t), intent(inout) :: this
1013
1014 if (this%is_malloced) call aero_rep_free_update_data(this%update_data)
1015
1016 end subroutine update_data_number_finalize
1017
1018!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1019
1020 !> Finalize an array of number update data objects
1022
1023 !> Update data objects to free
1024 type(aero_rep_update_data_single_particle_number_t), intent(inout) :: &
1025 this(:)
1026
1027 integer(kind=i_kind) :: i_aero
1028
1029 do i_aero = 1, size(this)
1030 call update_data_number_finalize(this(i_aero))
1031 end do
1032
1034
1035!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1036
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.
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