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)
79 type(
string_t),
allocatable,
private :: unique_names_(:)
81 type(
string_t),
allocatable,
private :: layer_names_(:)
83 integer(kind=i_kind) :: state_id_start = -99999
152 integer :: first_ = -9999
153 integer :: second_ = -9999
158 procedure :: constructor
166 logical :: is_malloced = .false.
168 integer(kind=i_kind) :: aero_rep_unique_id = 0
189 result(update_data)
bind (c)
192 type(c_ptr) :: update_data
197 update_data, aero_rep_unique_id, particle_id, number_conc) &
201 type(c_ptr),
value :: update_data
203 integer(kind=c_int),
value :: aero_rep_unique_id
205 integer(kind=c_int),
value :: particle_id
207 real(kind=c_double),
value :: number_conc
214 type(c_ptr),
value,
intent(in) :: update_data
249 integer(kind=i_kind),
intent(in) :: spec_state_id
252 type(
string_t),
allocatable :: layer_names_unordered(:)
254 type(
string_t),
allocatable :: cover_names_unordered(:)
257 integer(kind=i_kind),
allocatable :: ordered_layer_id(:)
259 character(len=:),
allocatable :: key_name, layer_name, layer_covers, &
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, &
270 key_name =
"maximum computational particles"
272 this%property_set%get_int(key_name, num_particles), &
273 "Missing maximum number of computational particles")
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//
"'")
286 allocate(phases(layers%size()))
287 allocate(cover_names_unordered(layers%size()))
288 allocate(layer_names_unordered(layers%size()))
293 call layers%iter_reset()
294 do i_layer = 1, layers%size()
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 '"// &
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
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 '"// &
317 cover_names_unordered(i_layer)%string = layer_covers
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 '"// &
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 '"// &
336 num_phases = num_phases + phases(i_layer)%val_%size()
338 call layers%iter_next()
343 cover_names_unordered)
346 allocate(this%layer_names_(
size(ordered_layer_id)))
347 this%layer_names_(:) = layer_names_unordered(ordered_layer_id(:))
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)
358 this%num_env_params = num_env_param_per_particle_ * num_particles
361 num_layers_ = layers%size()
362 max_particles_ = num_particles
367 allocate(this%aero_phase(num_phases * num_particles))
368 allocate(this%aero_phase_is_at_surface(num_phases * num_particles))
370 do i_layer = 1,
size(ordered_layer_id)
371 j_layer = ordered_layer_id(i_layer)
374 layer_phase_start_(i_layer) = curr_phase
375 layer_phase_end_(i_layer) = curr_phase + phases(j_layer)%val_%size() - 1
377 curr_phase = curr_phase + phases(j_layer)%val_%size()
383 do i_layer = 1,
size(ordered_layer_id)
384 j_layer = ordered_layer_id(i_layer)
386 call phases(j_layer)%val_%iter_reset()
387 do i_phase = 1, phases(j_layer)%val_%size()
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 '"// &
398 do j_phase = 1,
size(aero_phase_set)
399 if (aero_phase_set(j_phase)%val%name() .eq. phase_name)
then
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) = &
408 this%aero_phase_is_at_surface(i_particle*num_phases + curr_phase) = &
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
419 call assert(373124707, found)
421 call phases(j_layer)%val_%iter_next()
430 this%unique_names_ = this%unique_names( )
458 integer(kind=i_kind) :: state_size
462 state_size = max_particles_ * particle_state_size_
476 integer(kind=i_kind) :: state_size
480 state_size = particle_state_size_
503 character(len=*),
optional,
intent(in) :: phase_name
505 integer(kind=i_kind),
optional,
intent(in) :: tracer_type
507 character(len=*),
optional,
intent(in) ::
spec_name
509 logical,
optional,
intent(in) :: phase_is_at_surface
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(:)
516 if (.not.
present(phase_name) .and. &
517 .not.
present(tracer_type) .and. &
519 allocated(this%unique_names_))
then
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
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
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)
538 if (
spec_name .ne. spec_names(j_spec)%string) cycle
540 if (
present(tracer_type))
then
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
546 num_spec = num_spec + 1
548 deallocate(spec_names)
550 num_spec = num_spec + this%aero_phase(i_phase)%val%size()
555 num_spec = num_spec / max_particles_
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
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
567 spec_names = this%aero_phase(i_phase)%val%get_species_names()
568 do j_spec = 1, this%aero_phase(i_phase)%val%size()
570 if (
spec_name .ne. spec_names(j_spec)%string) cycle
572 if (
present(tracer_type))
then
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
578 do i_particle = 1, max_particles_
579 unique_names((i_particle-1)*num_spec+i_spec)%string =
'P'// &
581 this%layer_names_(i_layer)%string//
"."// &
582 this%aero_phase(i_phase)%val%name()//
"."// &
583 spec_names(j_spec)%string
587 deallocate(spec_names)
611 integer(kind=i_kind) :: spec_id
615 character(len=*),
intent(in) :: unique_name
617 integer(kind=i_kind) :: i_spec
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
624 spec_id = spec_id + 1
626 call die_msg( 449087541,
"Cannot find species '"//unique_name//
"'" )
637 character(len=:),
allocatable ::
spec_name
641 character(len=*),
intent(in) :: unique_name
643 type(string_t) :: l_unique_name
644 type(string_t),
allocatable :: substrs(:)
646 l_unique_name%string = unique_name
647 substrs = l_unique_name%split(
".")
648 call assert(407537518,
size( substrs ) .eq. 4 )
664 character(len=*),
intent(in) :: phase_name
666 logical,
intent(in),
optional :: is_at_surface
668 integer(kind=i_kind) :: i_phase, i_layer, phase_index
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
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
705 integer(kind=i_kind),
intent(in) :: phase_id
707 integer(kind=i_kind) :: i_phase
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 ) ) ) )
715 do i_phase = 1, total_num_phases_
717 this%aero_phase(i_phase)%val%num_jac_elem()
742 integer,
optional,
intent(in) :: layer
744 if (
present(layer))
then
745 num_phases = layer_phase_end_(layer) - layer_phase_start_(layer) + 1
747 num_phases = layer_phase_end_(num_layers_) - layer_phase_start_(1) + 1
762 integer,
optional,
intent(in) :: layer
764 integer,
optional,
intent(in) :: phase
766 if (
present(layer) .and.
present(phase))
then
767 if (layer .eq. num_layers_ .and. phase .eq. num_phases_(layer))
then
769 phase_state_id_(layer, phase)
770 else if (phase .eq. num_phases_(layer))
then
772 phase_state_id_(layer, phase)
775 phase_state_id_(layer, phase)
777 else if (
present(layer))
then
778 if (layer .eq. num_layers_)
then
780 phase_state_id_(layer, 1)
783 phase_state_id_(layer, 1)
785 else if (
present(phase))
then
786 call die_msg(917793122,
"Must specify layer if including phase is "// &
787 "state size request")
800 phase_name_second)
result(index_pairs)
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(:)
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
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_))
820 phase_id_first_ = -9999
821 phase_id_second_ = -9999
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
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
841 i_instance = i_instance + 1
846 allocate(temp_index_pairs(i_instance))
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
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
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_
888 if (
allocated(this%rep_name))
deallocate(this%rep_name)
889 if (
allocated(this%aero_phase))
then
891 call this%aero_phase(:)%dereference()
892 deallocate(this%aero_phase)
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)
912 integer(kind=i_kind) :: i_aero
914 do i_aero = 1,
size(this_array)
933 integer(kind=i_kind),
intent(in) :: aero_rep_type
936 if (aero_rep_id_.eq.-1)
then
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.
956 type(string_t),
intent(in) :: layer_names_unordered(:)
958 type(string_t),
intent(in) :: cover_names_unordered(:)
962 integer(kind=i_kind) :: i_layer, j_layer, i_cover
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)//
"'")
979 do i_layer = 1,
size(layer_names_unordered)
980 if (cover_names_unordered(i_layer)%string ==
"none")
then
987 do i_layer = 1,
size(layer_names_unordered)
989 .eq. cover_names_unordered(i_layer)%string)
then
1008 integer(kind=i_kind),
intent(in) :: particle_id
1010 real(kind=dp),
intent(in) :: number_conc
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, &
1033 integer,
intent(in) :: comm
1050 character,
intent(inout) :: buffer(:)
1052 integer,
intent(inout) :: pos
1054 integer,
intent(in) :: comm
1057 integer :: prev_position
1062 this%maximum_computational_particles, comm)
1064 call assert(411585487, &
1065 pos - prev_position <= this%pack_size(comm))
1078 character,
intent(inout) :: buffer(:)
1080 integer,
intent(inout) :: pos
1082 integer,
intent(in) :: comm
1085 integer :: prev_position
1090 this%maximum_computational_particles, comm)
1092 call assert(351557153, &
1093 pos - prev_position <= this%pack_size(comm))
1120 integer(kind=i_kind) :: i_aero
1122 do i_aero = 1,
size(this)
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...
Free an update data object.
Interface to c aerosol representation functions.
Set a new particle number concentration.
Interface for to_string functions.
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.
The chem_spec_data_t structure and associated subroutines.
Wrapper functions for MPI.
subroutine camp_mpi_pack_logical(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
subroutine camp_mpi_unpack_integer(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
subroutine camp_mpi_unpack_logical(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
integer function camp_mpi_pack_size_logical(val, comm)
Determines the number of bytes required to pack the given value.
subroutine camp_mpi_pack_integer(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
integer function camp_mpi_pack_size_integer(val, comm)
Determines the number of bytes required to pack the given value.
The property_t structure and associated subroutines.
Random number generators.
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.
Common utility subroutines.
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
subroutine die_msg(code, error_msg)
Error immediately.
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
character(len=camp_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Pointer type for building arrays.
Abstract aerosol representation data type.
Single particle aerosol representation.
Single particle update number concentration object.
String type for building arrays of string of various size.