PartMC
2.2.1
|
00001 ! Copyright (C) 2005-2012 Nicole Riemer and Matthew West 00002 ! Licensed under the GNU General Public License version 2 or (at your 00003 ! option) any later version. See the file COPYING for details. 00004 00005 !> \file 00006 !> The pmc_aero_state module. 00007 00008 !> The aero_state_t structure and assocated subroutines. 00009 module pmc_aero_state 00010 00011 use pmc_aero_particle_array 00012 use pmc_aero_sorted 00013 use pmc_integer_varray 00014 use pmc_bin_grid 00015 use pmc_aero_data 00016 use pmc_aero_particle 00017 use pmc_aero_dist 00018 use pmc_util 00019 use pmc_rand 00020 use pmc_aero_binned 00021 use pmc_mpi 00022 use pmc_spec_file 00023 use pmc_aero_info 00024 use pmc_aero_info_array 00025 use pmc_aero_weight 00026 #ifdef PMC_USE_MPI 00027 use mpi 00028 #endif 00029 00030 !> MPI tag for mixing particles between processes. 00031 integer, parameter :: AERO_STATE_TAG_MIX = 4987 00032 !> MPI tag for gathering between processes. 00033 integer, parameter :: AERO_STATE_TAG_GATHER = 4988 00034 !> MPI tag for scattering between processes. 00035 integer, parameter :: AERO_STATE_TAG_SCATTER = 4989 00036 00037 !> Single flat weighting scheme. 00038 integer, parameter :: AERO_STATE_WEIGHT_NONE = 1 00039 !> Single flat weighting scheme. 00040 integer, parameter :: AERO_STATE_WEIGHT_FLAT = 2 00041 !> Power-law weighting scheme. 00042 integer, parameter :: AERO_STATE_WEIGHT_POWER = 3 00043 !> Coupled number/mass weighting scheme. 00044 integer, parameter :: AERO_STATE_WEIGHT_NUMMASS = 4 00045 00046 !> The current collection of aerosol particles. 00047 !! 00048 !! The particles in \c aero_state_t are stored in a single flat 00049 !! array (the \c apa data member), with a sorting into bins possibly 00050 !! stored in the \c aero_sorted data member. 00051 !! 00052 !! Every time we remove particles we keep track of the particle ID 00053 !! and the action performed in the aero_info_array_t structure. This 00054 !! is typically cleared each time we output data to disk. 00055 type aero_state_t 00056 type(aero_particle_array_t) :: apa 00057 type(aero_sorted_t) :: aero_sorted 00058 logical :: valid_sort 00059 !> Weighting functions. 00060 type(aero_weight_t), allocatable :: aero_weight(:) 00061 !> Ideal number of computational particles. 00062 real(kind=dp) :: n_part_ideal 00063 !> Information on removed particles. 00064 type(aero_info_array_t) :: aero_info_array 00065 end type aero_state_t 00066 00067 contains 00068 00069 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00070 00071 !> Allocates aerosol arrays. 00072 subroutine aero_state_allocate(aero_state) 00073 00074 !> Aerosol to initialize. 00075 type(aero_state_t), intent(out) :: aero_state 00076 00077 call aero_particle_array_allocate(aero_state%apa) 00078 call aero_sorted_allocate(aero_state%aero_sorted) 00079 aero_state%valid_sort = .false. 00080 allocate(aero_state%aero_weight(0)) 00081 aero_state%n_part_ideal = 0d0 00082 call aero_info_array_allocate(aero_state%aero_info_array) 00083 00084 end subroutine aero_state_allocate 00085 00086 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00087 00088 !> Deallocates a previously allocated aerosol. 00089 subroutine aero_state_deallocate(aero_state) 00090 00091 !> Aerosol to deallocate. 00092 type(aero_state_t), intent(inout) :: aero_state 00093 00094 call aero_particle_array_deallocate(aero_state%apa) 00095 call aero_sorted_deallocate(aero_state%aero_sorted) 00096 aero_state%valid_sort = .false. 00097 call aero_weight_deallocate(aero_state%aero_weight) 00098 deallocate(aero_state%aero_weight) 00099 call aero_info_array_deallocate(aero_state%aero_info_array) 00100 00101 end subroutine aero_state_deallocate 00102 00103 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00104 00105 !> Resets an \c aero_state to an empty state. 00106 subroutine aero_state_reset(aero_state) 00107 00108 !> Aerosol to reset. 00109 type(aero_state_t), intent(inout) :: aero_state 00110 00111 call aero_state_deallocate(aero_state) 00112 call aero_state_allocate(aero_state) 00113 00114 end subroutine aero_state_reset 00115 00116 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00117 00118 !> Copies aerosol to a destination that has already had 00119 !> aero_state_allocate() called on it. 00120 subroutine aero_state_copy(aero_state_from, aero_state_to) 00121 00122 !> Reference aerosol. 00123 type(aero_state_t), intent(in) :: aero_state_from 00124 !> Already allocated. 00125 type(aero_state_t), intent(inout) :: aero_state_to 00126 00127 call aero_particle_array_copy(aero_state_from%apa, aero_state_to%apa) 00128 aero_state_to%valid_sort = .false. 00129 call aero_state_copy_weight(aero_state_from, aero_state_to) 00130 aero_state_to%n_part_ideal = aero_state_from%n_part_ideal 00131 call aero_info_array_copy(aero_state_from%aero_info_array, & 00132 aero_state_to%aero_info_array) 00133 00134 end subroutine aero_state_copy 00135 00136 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00137 00138 !> Copies weighting information for an \c aero_state. 00139 subroutine aero_state_copy_weight(aero_state_from, aero_state_to) 00140 00141 !> Reference aerosol. 00142 type(aero_state_t), intent(in) :: aero_state_from 00143 !> Already allocated. 00144 type(aero_state_t), intent(inout) :: aero_state_to 00145 00146 call aero_weight_array_copy(aero_state_from%aero_weight, & 00147 aero_state_to%aero_weight) 00148 00149 end subroutine aero_state_copy_weight 00150 00151 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00152 00153 !> Sets the weighting functions for an \c aero_state. 00154 subroutine aero_state_set_weight(aero_state, weight_type, exponent) 00155 00156 !> Aerosol to set the weights on. 00157 type(aero_state_t), intent(inout) :: aero_state 00158 !> Type of weighting scheme to use. 00159 integer, intent(in) :: weight_type 00160 !> Exponent for power-law weighting (only used if \c weight_type 00161 !> is \c AERO_STATE_WEIGHT_POWER). 00162 real(kind=dp), intent(in), optional :: exponent 00163 00164 aero_state%valid_sort = .false. 00165 call aero_weight_deallocate(aero_state%aero_weight) 00166 deallocate(aero_state%aero_weight) 00167 if (weight_type == AERO_STATE_WEIGHT_NONE) then 00168 allocate(aero_state%aero_weight(0)) 00169 elseif (weight_type == AERO_STATE_WEIGHT_FLAT) then 00170 allocate(aero_state%aero_weight(1)) 00171 call aero_weight_allocate(aero_state%aero_weight(1)) 00172 aero_state%aero_weight(1)%type = AERO_WEIGHT_TYPE_NONE 00173 aero_state%aero_weight(1)%ref_radius = 1d0 00174 aero_state%aero_weight(1)%exponent = 0d0 00175 elseif (weight_type == AERO_STATE_WEIGHT_POWER) then 00176 call assert_msg(656670336, present(exponent), & 00177 "exponent parameter required for AERO_STATE_WEIGHT_POWER") 00178 allocate(aero_state%aero_weight(1)) 00179 call aero_weight_allocate(aero_state%aero_weight(1)) 00180 aero_state%aero_weight(1)%type = AERO_WEIGHT_TYPE_POWER 00181 aero_state%aero_weight(1)%ref_radius = 1d0 00182 aero_state%aero_weight(1)%exponent = exponent 00183 elseif (weight_type == AERO_STATE_WEIGHT_NUMMASS) then 00184 allocate(aero_state%aero_weight(2)) 00185 call aero_weight_allocate(aero_state%aero_weight(1)) 00186 aero_state%aero_weight(1)%type = AERO_WEIGHT_TYPE_NONE 00187 aero_state%aero_weight(1)%ref_radius = 1d0 00188 aero_state%aero_weight(1)%exponent = 0d0 00189 call aero_weight_allocate(aero_state%aero_weight(2)) 00190 aero_state%aero_weight(2)%type = AERO_WEIGHT_TYPE_POWER 00191 aero_state%aero_weight(2)%ref_radius = 1d0 00192 aero_state%aero_weight(2)%exponent = -3d0 00193 else 00194 call die_msg(969076992, "unknown weight_type: " & 00195 // trim(integer_to_string(weight_type))) 00196 end if 00197 00198 end subroutine aero_state_set_weight 00199 00200 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00201 00202 !> Returns the total number of particles in an aerosol distribution. 00203 integer function aero_state_total_particles(aero_state, i_group) 00204 00205 !> Aerosol state. 00206 type(aero_state_t), intent(in) :: aero_state 00207 !> Weight group. 00208 integer, optional, intent(in) :: i_group 00209 00210 integer :: i_part 00211 00212 if (present(i_group)) then 00213 if (aero_state%valid_sort) then 00214 aero_state_total_particles & 00215 = aero_state%aero_sorted%weight%inverse(i_group)%n_entry 00216 else 00217 ! FIXME: should we just sort? 00218 aero_state_total_particles = 0 00219 do i_part = 1,aero_state%apa%n_part 00220 if (aero_state%apa%particle(i_part)%weight_group == i_group) then 00221 aero_state_total_particles = aero_state_total_particles + 1 00222 end if 00223 end do 00224 end if 00225 else 00226 aero_state_total_particles = aero_state%apa%n_part 00227 end if 00228 00229 end function aero_state_total_particles 00230 00231 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00232 00233 !> Returns the total number of particles across all processes. 00234 integer function aero_state_total_particles_all_procs(aero_state, i_group) 00235 00236 !> Aerosol state. 00237 type(aero_state_t), intent(in) :: aero_state 00238 !> Weight group. 00239 integer, optional, intent(in) :: i_group 00240 00241 call pmc_mpi_allreduce_sum_integer(& 00242 aero_state_total_particles(aero_state, i_group), & 00243 aero_state_total_particles_all_procs) 00244 00245 end function aero_state_total_particles_all_procs 00246 00247 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00248 00249 !> Resets an aero_state to have zero particles per bin. This must 00250 !> already have had aero_state_allocate() called on it. This 00251 !> function can be called more than once on the same state. 00252 subroutine aero_state_zero(aero_state) 00253 00254 !> State to zero. 00255 type(aero_state_t), intent(inout) :: aero_state 00256 00257 integer :: i, n_bin 00258 00259 call aero_particle_array_zero(aero_state%apa) 00260 aero_state%valid_sort = .false. 00261 call aero_info_array_zero(aero_state%aero_info_array) 00262 00263 end subroutine aero_state_zero 00264 00265 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00266 00267 !> Add the given particle. 00268 subroutine aero_state_add_particle(aero_state, aero_particle, allow_resort) 00269 00270 !> Aerosol state. 00271 type(aero_state_t), intent(inout) :: aero_state 00272 !> Particle to add. 00273 type(aero_particle_t), intent(in) :: aero_particle 00274 !> Whether to allow a resort due to the add. 00275 logical, optional, intent(in) :: allow_resort 00276 00277 if (aero_state%valid_sort) then 00278 call aero_sorted_add_particle(aero_state%aero_sorted, aero_state%apa, & 00279 aero_particle, size(aero_state%aero_weight), allow_resort) 00280 else 00281 call aero_particle_array_add_particle(aero_state%apa, aero_particle) 00282 end if 00283 00284 end subroutine aero_state_add_particle 00285 00286 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00287 00288 !> Remove the given particle without recording it. 00289 subroutine aero_state_remove_particle_no_info(aero_state, i_part) 00290 00291 !> Aerosol state. 00292 type(aero_state_t), intent(inout) :: aero_state 00293 !> Index of particle to remove. 00294 integer, intent(in) :: i_part 00295 00296 if (aero_state%valid_sort) then 00297 call aero_sorted_remove_particle(aero_state%aero_sorted, & 00298 aero_state%apa, i_part) 00299 else 00300 call aero_particle_array_remove_particle(aero_state%apa, i_part) 00301 end if 00302 00303 end subroutine aero_state_remove_particle_no_info 00304 00305 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00306 00307 !> Remove the given particle and record the removal. 00308 subroutine aero_state_remove_particle_with_info(aero_state, i_part, & 00309 aero_info) 00310 00311 !> Aerosol state. 00312 type(aero_state_t), intent(inout) :: aero_state 00313 !> Index of particle to remove. 00314 integer, intent(in) :: i_part 00315 !> Removal info. 00316 type(aero_info_t), intent(in) :: aero_info 00317 00318 call aero_state_remove_particle_no_info(aero_state, i_part) 00319 call aero_info_array_add_aero_info(aero_state%aero_info_array, & 00320 aero_info) 00321 00322 end subroutine aero_state_remove_particle_with_info 00323 00324 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00325 00326 !> Remove the given particle and possibly record the removal. 00327 subroutine aero_state_remove_particle(aero_state, i_part, & 00328 record_removal, aero_info) 00329 00330 !> Aerosol state. 00331 type(aero_state_t), intent(inout) :: aero_state 00332 !> Index of particle to remove. 00333 integer, intent(in) :: i_part 00334 !> Whether to record the removal in the aero_info_array. 00335 logical, intent(in) :: record_removal 00336 !> Removal info. 00337 type(aero_info_t), intent(in) :: aero_info 00338 00339 if (record_removal) then 00340 call aero_state_remove_particle_with_info(aero_state, i_part, & 00341 aero_info) 00342 else 00343 call aero_state_remove_particle_no_info(aero_state, i_part) 00344 end if 00345 00346 end subroutine aero_state_remove_particle 00347 00348 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00349 00350 !> Remove a randomly chosen particle from the given bin and return 00351 !> it. 00352 subroutine aero_state_remove_rand_particle_from_bin(aero_state, & 00353 i_bin, aero_particle) 00354 00355 !> Aerosol state. 00356 type(aero_state_t), intent(inout) :: aero_state 00357 !> Bin number to remove particle from. 00358 integer, intent(in) :: i_bin 00359 !> Removed particle. 00360 type(aero_particle_t), intent(inout) :: aero_particle 00361 00362 integer :: i_entry, i_part 00363 00364 call assert(742996300, aero_state%valid_sort) 00365 call assert(392182617, & 00366 aero_state%aero_sorted%size%inverse(i_bin)%n_entry > 0) 00367 i_entry = pmc_rand_int(aero_state%aero_sorted%size%inverse(i_bin)%n_entry) 00368 i_part = aero_state%aero_sorted%size%inverse(i_bin)%entry(i_entry) 00369 call aero_particle_copy(aero_state%apa%particle(i_part), aero_particle) 00370 call aero_state_remove_particle_no_info(aero_state, i_part) 00371 00372 end subroutine aero_state_remove_rand_particle_from_bin 00373 00374 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00375 00376 !> Add copies or remove a particle, with a given mean number of 00377 !> resulting particles. 00378 !! 00379 !! The particle number \c i_part is either removed, or zero or more 00380 !! copies are added, with a random number of copies with the given 00381 !! mean \c n_part_mean. The final number of particles is either 00382 !! <tt>floor(n_part_mean)</tt> or <tt>ceiling(n_part_mean)</tt>, 00383 !! chosen randomly so the mean is \c n_part_mean. 00384 subroutine aero_state_dup_particle(aero_state, i_part, n_part_mean, & 00385 random_weight_group) 00386 00387 !> Aerosol state. 00388 type(aero_state_t), intent(inout) :: aero_state 00389 !> Particle number. 00390 integer, intent(in) :: i_part 00391 !> Mean number of resulting particles. 00392 real(kind=dp), intent(in) :: n_part_mean 00393 !> Whether particle copies should be placed in a randomly chosen 00394 !> weight group. 00395 logical, optional, intent(in) :: random_weight_group 00396 00397 integer :: n_copies, i_dup, new_group 00398 type(aero_particle_t), pointer :: aero_particle 00399 type(aero_particle_t) :: new_aero_particle 00400 type(aero_info_t) :: aero_info 00401 00402 aero_particle => aero_state%apa%particle(i_part) 00403 n_copies = prob_round(n_part_mean) 00404 if (n_copies == 0) then 00405 call aero_info_allocate(aero_info) 00406 aero_info%id = aero_particle%id 00407 aero_info%action = AERO_INFO_WEIGHT 00408 aero_info%other_id = 0 00409 call aero_state_remove_particle_with_info(aero_state, & 00410 i_part, aero_info) 00411 call aero_info_deallocate(aero_info) 00412 elseif (n_copies > 1) then 00413 call aero_particle_allocate(new_aero_particle) 00414 do i_dup = 1,(n_copies - 1) 00415 call aero_particle_copy(aero_particle, new_aero_particle) 00416 call aero_particle_new_id(new_aero_particle) 00417 if (present(random_weight_group)) then 00418 if (random_weight_group) then 00419 new_group & 00420 = aero_weight_array_rand_group(aero_state%aero_weight, & 00421 aero_particle_radius(aero_particle)) 00422 call aero_particle_set_group(new_aero_particle, new_group) 00423 end if 00424 end if 00425 call aero_state_add_particle(aero_state, new_aero_particle) 00426 ! re-get the particle pointer, which may have 00427 ! changed due to reallocations caused by adding 00428 aero_particle => aero_state%apa%particle(i_part) 00429 end do 00430 call aero_particle_deallocate(new_aero_particle) 00431 end if 00432 00433 end subroutine aero_state_dup_particle 00434 00435 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00436 00437 !> The number concentration of a single particle (m^{-3}). 00438 real(kind=dp) function aero_state_particle_num_conc(aero_state, & 00439 aero_particle) 00440 00441 !> Aerosol state containing the particle. 00442 type(aero_state_t), intent(in) :: aero_state 00443 !> Aerosol particle. 00444 type(aero_particle_t), intent(in) :: aero_particle 00445 00446 aero_state_particle_num_conc & 00447 = aero_weight_array_num_conc(aero_state%aero_weight, aero_particle) 00448 00449 end function aero_state_particle_num_conc 00450 00451 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00452 00453 !> Save the correct number concentrations for later use by 00454 !> aero_state_reweight(). 00455 subroutine aero_state_num_conc_for_reweight(aero_state, reweight_num_conc) 00456 00457 !> Aerosol state. 00458 type(aero_state_t), intent(in) :: aero_state 00459 !> Number concentrations for later use by aero_state_reweight(). 00460 real(kind=dp), intent(out) :: reweight_num_conc(aero_state%apa%n_part) 00461 00462 integer :: i_part 00463 00464 do i_part = 1,aero_state%apa%n_part 00465 reweight_num_conc(i_part) & 00466 = aero_weight_array_single_num_conc(aero_state%aero_weight, & 00467 aero_state%apa%particle(i_part)) 00468 end do 00469 00470 end subroutine aero_state_num_conc_for_reweight 00471 00472 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00473 00474 !> Reweight all particles after their constituent volumes have been 00475 !> altered. 00476 !! 00477 !! The pattern for use should be like: 00478 !! <pre> 00479 !! call aero_state_num_conc_for_reweight(aero_state, reweight_num_conc) 00480 !! ... alter particle species volumes in aero_state ... 00481 !! call aero_state_reweight(aero_state, reweight_num_conc) 00482 !! </pre> 00483 subroutine aero_state_reweight(aero_state, reweight_num_conc) 00484 00485 !> Aerosol state. 00486 type(aero_state_t), intent(inout) :: aero_state 00487 !> Number concentrations previously computed by 00488 !> aero_state_num_conc_for_reweight(). 00489 real(kind=dp), intent(in) :: reweight_num_conc(aero_state%apa%n_part) 00490 00491 integer :: i_part, i_group 00492 real(kind=dp) :: n_part_old(size(aero_state%aero_weight)) 00493 real(kind=dp) :: n_part_new(size(aero_state%aero_weight)) 00494 real(kind=dp) :: old_num_conc, new_num_conc, n_part_mean 00495 type(aero_particle_t), pointer :: aero_particle 00496 00497 ! find average number of new particles in each weight group, if 00498 ! comp_vol is not changed 00499 n_part_old = 0d0 00500 n_part_new = 0d0 00501 do i_part = 1,aero_state%apa%n_part 00502 aero_particle => aero_state%apa%particle(i_part) 00503 old_num_conc = reweight_num_conc(i_part) 00504 new_num_conc & 00505 = aero_weight_array_single_num_conc(aero_state%aero_weight, & 00506 aero_particle) 00507 n_part_mean = old_num_conc / new_num_conc 00508 n_part_new(aero_particle%weight_group) & 00509 = n_part_new(aero_particle%weight_group) & 00510 + n_part_mean 00511 n_part_old(aero_particle%weight_group) & 00512 = n_part_old(aero_particle%weight_group) + 1d0 00513 end do 00514 00515 ! alter comp_vol to leave the number of computational particles 00516 ! per weight bin unchanged 00517 do i_group = 1,size(aero_state%aero_weight) 00518 if (n_part_old(i_group) == 0d0) cycle 00519 aero_state%aero_weight(i_group)%comp_vol & 00520 = aero_state%aero_weight(i_group)%comp_vol & 00521 * n_part_old(i_group) / n_part_new(i_group) 00522 end do 00523 00524 ! work backwards so any additions and removals will only affect 00525 ! particles that we've already dealt with 00526 do i_part = aero_state%apa%n_part,1,-1 00527 aero_particle => aero_state%apa%particle(i_part) 00528 old_num_conc = reweight_num_conc(i_part) 00529 new_num_conc & 00530 = aero_weight_array_single_num_conc(aero_state%aero_weight, & 00531 aero_particle) 00532 n_part_mean = old_num_conc / new_num_conc 00533 call aero_state_dup_particle(aero_state, i_part, n_part_mean) 00534 end do 00535 00536 end subroutine aero_state_reweight 00537 00538 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00539 00540 !> <tt>aero_state += aero_state_delta</tt>, with adding the computational 00541 !> volumes, so the new concentration is the (volume-weighted) 00542 !> average of the two concentration. 00543 subroutine aero_state_add(aero_state, aero_state_delta) 00544 00545 !> Aerosol state. 00546 type(aero_state_t), intent(inout) :: aero_state 00547 !> Increment. 00548 type(aero_state_t), intent(in) :: aero_state_delta 00549 00550 call aero_state_add_particles(aero_state, aero_state_delta) 00551 call aero_weight_add_comp_vol(aero_state%aero_weight, & 00552 aero_state_delta%aero_weight) 00553 00554 end subroutine aero_state_add 00555 00556 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00557 00558 !> <tt>aero_state += aero_state_delta</tt>, with the computational 00559 !> volume of \c aero_state left unchanged, so the new concentration is the 00560 !> sum of the two concentrations, computed with \c aero_state%comp_vol. 00561 subroutine aero_state_add_particles(aero_state, aero_state_delta) 00562 00563 !> Aerosol state. 00564 type(aero_state_t), intent(inout) :: aero_state 00565 !> Increment. 00566 type(aero_state_t), intent(in) :: aero_state_delta 00567 00568 integer :: i_part, i_bin 00569 00570 do i_part = 1,aero_state_delta%apa%n_part 00571 call aero_state_add_particle(aero_state, & 00572 aero_state_delta%apa%particle(i_part)) 00573 end do 00574 call aero_info_array_add(aero_state%aero_info_array, & 00575 aero_state_delta%aero_info_array) 00576 00577 end subroutine aero_state_add_particles 00578 00579 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00580 00581 !> Generates a Poisson sample of an \c aero_dist, adding to \c 00582 !> aero_state. The sampled amount is <tt>sample_prop * 00583 !> aero_state%comp_vol</tt>. 00584 subroutine aero_state_add_aero_dist_sample(aero_state, aero_data, & 00585 aero_dist, sample_prop, create_time, n_part_add) 00586 00587 !> Aero state to add to. 00588 type(aero_state_t), intent(inout) :: aero_state 00589 !> Aero data values. 00590 type(aero_data_t), intent(in) :: aero_data 00591 !> Distribution to sample. 00592 type(aero_dist_t), intent(in) :: aero_dist 00593 !> Volume fraction to sample (1). 00594 real(kind=dp), intent(in) :: sample_prop 00595 !> Creation time for new particles (s). 00596 real(kind=dp), intent(in) :: create_time 00597 !> Number of particles added. 00598 integer, intent(out), optional :: n_part_add 00599 00600 real(kind=dp) :: n_samp_avg, radius, vol, mean_n_part, n_part_new 00601 real(kind=dp) :: comp_vol_ratio, n_part_ideal_local_group 00602 integer :: n_samp, i_mode, i_samp, i_group, n_group, global_n_part 00603 type(aero_mode_t), pointer :: aero_mode 00604 type(aero_particle_t) :: aero_particle 00605 00606 call aero_particle_allocate_size(aero_particle, aero_data%n_spec, & 00607 aero_data%n_source) 00608 00609 n_group = size(aero_state%aero_weight) 00610 if (present(n_part_add)) then 00611 n_part_add = 0 00612 end if 00613 do i_group = 1,n_group 00614 global_n_part & 00615 = aero_state_total_particles_all_procs(aero_state, i_group) 00616 mean_n_part = real(global_n_part, kind=dp) 00617 / real(pmc_mpi_size(), kind=dp) 00618 if (aero_state%aero_weight(i_group)%comp_vol == 0d0) then 00619 ! FIXME: assert that n_part in this weight group is zero 00620 aero_state%aero_weight(i_group)%comp_vol = 1d0 00621 end if 00622 n_samp_avg = sample_prop & 00623 * aero_dist_number(aero_dist, aero_state%aero_weight(i_group)) 00624 n_part_new = mean_n_part + n_samp_avg 00625 if (n_part_new == 0d0) cycle 00626 n_part_ideal_local_group = aero_state%n_part_ideal & 00627 / real(pmc_mpi_size(), kind=dp) / real(n_group, kind=dp) 00628 if ((n_part_new < n_part_ideal_local_group / 2d0) & 00629 .or. (n_part_new > n_part_ideal_local_group * 2d0)) & 00630 then 00631 comp_vol_ratio = n_part_ideal_local_group / n_part_new 00632 call aero_state_scale_comp_vol(aero_state, i_group, comp_vol_ratio) 00633 end if 00634 do i_mode = 1,aero_dist%n_mode 00635 aero_mode => aero_dist%mode(i_mode) 00636 n_samp_avg = sample_prop & 00637 * aero_mode_number(aero_mode, aero_state%aero_weight(i_group)) 00638 n_samp = rand_poisson(n_samp_avg) 00639 if (present(n_part_add)) then 00640 n_part_add = n_part_add + n_samp 00641 end if 00642 do i_samp = 1,n_samp 00643 call aero_particle_zero(aero_particle) 00644 call aero_mode_sample_radius(aero_mode, & 00645 aero_state%aero_weight(i_group), radius) 00646 vol = rad2vol(radius) 00647 call aero_particle_set_vols(aero_particle, & 00648 aero_mode%vol_frac * vol) 00649 call aero_particle_new_id(aero_particle) 00650 call aero_particle_set_group(aero_particle, i_group) 00651 call aero_particle_set_create_time(aero_particle, create_time) 00652 call aero_particle_set_source(aero_particle, aero_mode%source) 00653 call aero_state_add_particle(aero_state, aero_particle) 00654 end do 00655 end do 00656 end do 00657 call aero_particle_deallocate(aero_particle) 00658 00659 end subroutine aero_state_add_aero_dist_sample 00660 00661 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00662 00663 !> Choose a random particle from the aero_state. 00664 subroutine aero_state_rand_particle(aero_state, i_part) 00665 00666 !> Original state. 00667 type(aero_state_t), intent(in) :: aero_state 00668 !> Chosen random particle number. 00669 integer, intent(out) :: i_part 00670 00671 call assert(950725003, aero_state%apa%n_part > 0) 00672 i_part = pmc_rand_int(aero_state%apa%n_part) 00673 00674 end subroutine aero_state_rand_particle 00675 00676 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00677 00678 !> Generates a random sample by removing particles from 00679 !> aero_state_from and adding them to aero_state_to, which must 00680 !> be already allocated (and should have its comp_vol set). 00681 !! 00682 !! None of the computational volumes are altered by this sampling, 00683 !! making this the equivalent of aero_state_add_particles(). 00684 subroutine aero_state_sample_particles(aero_state_from, aero_state_to, & 00685 sample_prob, removal_action) 00686 00687 !> Original state. 00688 type(aero_state_t), intent(inout) :: aero_state_from 00689 !> Destination state. 00690 type(aero_state_t), intent(inout) :: aero_state_to 00691 !> Probability of sampling each particle. 00692 real(kind=dp), intent(in) :: sample_prob 00693 !> Action for removal (see pmc_aero_info module for action 00694 !> parameters). Set to AERO_INFO_NONE to not log removal. 00695 integer, intent(in) :: removal_action 00696 00697 integer :: n_transfer, i_transfer, i_part 00698 logical :: do_add, do_remove 00699 real(kind=dp) :: num_conc_from, num_conc_to 00700 type(aero_info_t) :: aero_info 00701 00702 call assert(721006962, (sample_prob >= 0d0) .and. (sample_prob <= 1d0)) 00703 call aero_state_reset(aero_state_to) 00704 call aero_state_copy_weight(aero_state_from, aero_state_to) 00705 n_transfer = rand_binomial(aero_state_total_particles(aero_state_from), & 00706 sample_prob) 00707 i_transfer = 0 00708 do while (i_transfer < n_transfer) 00709 if (aero_state_total_particles(aero_state_from) <= 0) exit 00710 call aero_state_rand_particle(aero_state_from, i_part) 00711 num_conc_from = aero_weight_array_num_conc( & 00712 aero_state_from%aero_weight, & 00713 aero_state_from%apa%particle(i_part)) 00714 num_conc_to = aero_weight_array_num_conc( & 00715 aero_state_to%aero_weight, & 00716 aero_state_from%apa%particle(i_part)) 00717 00718 if (num_conc_to == num_conc_from) then ! add and remove 00719 do_add = .true. 00720 do_remove = .true. 00721 elseif (num_conc_to < num_conc_from) then ! always add, maybe remove 00722 do_add = .true. 00723 do_remove = .false. 00724 if (pmc_random() < num_conc_to / num_conc_from) then 00725 do_remove = .true. 00726 end if 00727 else ! always remove, maybe add 00728 do_add = .false. 00729 if (pmc_random() < num_conc_from / num_conc_to) then 00730 do_add = .true. 00731 end if 00732 do_remove = .true. 00733 end if 00734 if (do_add) then 00735 call aero_state_add_particle(aero_state_to, & 00736 aero_state_from%apa%particle(i_part)) 00737 end if 00738 if (do_remove) then 00739 if (removal_action /= AERO_INFO_NONE) then 00740 call aero_info_allocate(aero_info) 00741 aero_info%id = aero_state_from%apa%particle(i_part)%id 00742 aero_info%action = removal_action 00743 call aero_state_remove_particle_with_info(aero_state_from, & 00744 i_part, aero_info) 00745 call aero_info_deallocate(aero_info) 00746 else 00747 call aero_state_remove_particle_no_info(aero_state_from, & 00748 i_part) 00749 end if 00750 i_transfer = i_transfer + 1 00751 end if 00752 end do 00753 00754 end subroutine aero_state_sample_particles 00755 00756 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00757 00758 !> Generates a random sample by removing particles from 00759 !> aero_state_from and adding them to aero_state_to, transfering 00760 !> computational volume as well. This is the equivalent of 00761 !> aero_state_add(). 00762 subroutine aero_state_sample(aero_state_from, aero_state_to, & 00763 sample_prob, removal_action) 00764 00765 !> Original state. 00766 type(aero_state_t), intent(inout) :: aero_state_from 00767 !> Destination state (previous contents will be lost). 00768 type(aero_state_t), intent(inout) :: aero_state_to 00769 !> Probability of sampling each particle. 00770 real(kind=dp), intent(in) :: sample_prob 00771 !> Action for removal (see pmc_aero_info module for action 00772 !> parameters). Set to AERO_INFO_NONE to not log removal. 00773 integer, intent(in) :: removal_action 00774 00775 integer :: n_transfer, i_transfer, i_part 00776 logical :: do_add, do_remove 00777 real(kind=dp) :: num_conc_from, num_conc_to 00778 type(aero_info_t) :: aero_info 00779 00780 call assert(393205561, (sample_prob >= 0d0) .and. (sample_prob <= 1d0)) 00781 call aero_state_reset(aero_state_to) 00782 call aero_state_copy_weight(aero_state_from, aero_state_to) 00783 call aero_weight_zero_comp_vol(aero_state_to%aero_weight) 00784 n_transfer = rand_binomial(aero_state_total_particles(aero_state_from), & 00785 sample_prob) 00786 do i_transfer = 1,n_transfer 00787 if (aero_state_total_particles(aero_state_from) <= 0) exit 00788 call aero_state_rand_particle(aero_state_from, i_part) 00789 00790 call aero_state_add_particle(aero_state_to, & 00791 aero_state_from%apa%particle(i_part)) 00792 if (removal_action /= AERO_INFO_NONE) then 00793 call aero_info_allocate(aero_info) 00794 aero_info%id = aero_state_from%apa%particle(i_part)%id 00795 aero_info%action = removal_action 00796 call aero_state_remove_particle_with_info(aero_state_from, & 00797 i_part, aero_info) 00798 call aero_info_deallocate(aero_info) 00799 else 00800 call aero_state_remove_particle_no_info(aero_state_from, & 00801 i_part) 00802 end if 00803 end do 00804 call aero_weight_transfer_comp_vol(aero_state_from%aero_weight, & 00805 aero_state_to%aero_weight, sample_prob) 00806 00807 end subroutine aero_state_sample 00808 00809 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00810 00811 !> Create the bin number and mass arrays from aero_state%v. 00812 subroutine aero_state_to_binned(bin_grid, aero_data, aero_state, & 00813 aero_binned) 00814 00815 !> Bin grid. 00816 type(bin_grid_t), intent(in) :: bin_grid 00817 !> Aerosol data. 00818 type(aero_data_t), intent(in) :: aero_data 00819 !> Aerosol state. 00820 type(aero_state_t), intent(in) :: aero_state 00821 !> Binned distributions. 00822 type(aero_binned_t), intent(inout) :: aero_binned 00823 00824 integer :: i_part, i_bin 00825 type(aero_particle_t), pointer :: aero_particle 00826 00827 aero_binned%num_conc = 0d0 00828 aero_binned%vol_conc = 0d0 00829 do i_part = 1,aero_state%apa%n_part 00830 aero_particle => aero_state%apa%particle(i_part) 00831 i_bin = bin_grid_particle_in_bin(bin_grid, & 00832 aero_particle_radius(aero_particle)) 00833 if ((i_bin < 1) .or. (i_bin > bin_grid%n_bin)) then 00834 call warn_msg(980232449, "particle ID " & 00835 // trim(integer_to_string(aero_particle%id)) & 00836 // " outside of bin_grid, discarding") 00837 else 00838 aero_binned%vol_conc(i_bin,:) = aero_binned%vol_conc(i_bin,:) & 00839 + aero_particle%vol & 00840 * aero_weight_array_num_conc(aero_state%aero_weight, & 00841 aero_particle) / bin_grid%log_width 00842 aero_binned%num_conc(i_bin) = aero_binned%num_conc(i_bin) & 00843 + aero_weight_array_num_conc(aero_state%aero_weight, & 00844 aero_particle) / bin_grid%log_width 00845 end if 00846 end do 00847 00848 end subroutine aero_state_to_binned 00849 00850 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00851 00852 !> Returns the diameters of all particles. The \c diameters array 00853 !> will be reallocated if necessary. 00854 subroutine aero_state_diameters(aero_state, diameters) 00855 00856 !> Aerosol state. 00857 type(aero_state_t), intent(in) :: aero_state 00858 !> Diameters array (m). 00859 real(kind=dp), intent(inout), allocatable :: diameters(:) 00860 00861 call ensure_real_array_size(diameters, aero_state%apa%n_part) 00862 diameters = aero_particle_diameter( & 00863 aero_state%apa%particle(1:aero_state%apa%n_part)) 00864 00865 end subroutine aero_state_diameters 00866 00867 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00868 00869 !> Returns the masses of all particles. The \c masses array 00870 !> will be reallocated if necessary. 00871 !! 00872 !! If \c include is specified then only those species are included 00873 !! in computing the masses. If \c exclude is specified then all 00874 !! species except those species are included. If both \c include and 00875 !! \c exclude arguments are specified then only those species in \c 00876 !! include but not in \c exclude are included. 00877 subroutine aero_state_masses(aero_state, aero_data, masses, include, exclude) 00878 00879 !> Aerosol state. 00880 type(aero_state_t), intent(in) :: aero_state 00881 !> Aerosol data. 00882 type(aero_data_t), intent(in) :: aero_data 00883 !> Masses array (kg). 00884 real(kind=dp), intent(inout), allocatable :: masses(:) 00885 !> Species names to include in the mass. 00886 character(len=*), optional :: include(:) 00887 !> Species names to exclude in the mass. 00888 character(len=*), optional :: exclude(:) 00889 00890 logical :: use_species(aero_data%n_spec) 00891 integer :: i_name, i_spec 00892 00893 call ensure_real_array_size(masses, aero_state%apa%n_part) 00894 if ((.not. present(include)) .and. (.not. present(exclude))) then 00895 masses = aero_particle_mass( & 00896 aero_state%apa%particle(1:aero_state%apa%n_part), aero_data) 00897 else 00898 if (present(include)) then 00899 use_species = .false. 00900 do i_name = 1, size(include) 00901 i_spec = aero_data_spec_by_name(aero_data, include(i_name)) 00902 call assert_msg(963163690, i_spec > 0, & 00903 "unknown species: " // trim(include(i_name))) 00904 use_species(i_spec) = .true. 00905 end do 00906 else 00907 use_species = .true. 00908 end if 00909 if (present(exclude)) then 00910 do i_name = 1, size(exclude) 00911 i_spec = aero_data_spec_by_name(aero_data, exclude(i_name)) 00912 call assert_msg(950847713, i_spec > 0, & 00913 "unknown species: " // trim(exclude(i_name))) 00914 use_species(i_spec) = .false. 00915 end do 00916 end if 00917 masses = 0d0 00918 do i_spec = 1,aero_data%n_spec 00919 if (use_species(i_spec)) then 00920 masses = masses + aero_particle_species_mass( & 00921 aero_state%apa%particle(1:aero_state%apa%n_part), & 00922 i_spec, aero_data) 00923 end if 00924 end do 00925 end if 00926 00927 end subroutine aero_state_masses 00928 00929 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00930 00931 !> Returns the number concentrations of all particles. The \c 00932 !> num_concs array will be reallocated if necessary. 00933 subroutine aero_state_num_concs(aero_state, num_concs) 00934 00935 !> Aerosol state. 00936 type(aero_state_t), intent(in) :: aero_state 00937 !> Number concentrations array (m^{-3}). 00938 real(kind=dp), intent(inout), allocatable :: num_concs(:) 00939 00940 integer :: i_part 00941 00942 call ensure_real_array_size(num_concs, aero_state%apa%n_part) 00943 do i_part = 1,aero_state%apa%n_part 00944 num_concs(i_part) = aero_state_particle_num_conc(aero_state, & 00945 aero_state%apa%particle(i_part)) 00946 end do 00947 00948 end subroutine aero_state_num_concs 00949 00950 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00951 00952 !> Returns the total number concentration. 00953 real(kind=dp) function aero_state_total_num_conc(aero_state) 00954 00955 !> Aerosol state. 00956 type(aero_state_t), intent(in) :: aero_state 00957 00958 integer :: i_part 00959 00960 aero_state_total_num_conc = 0d0 00961 do i_part = 1,aero_state%apa%n_part 00962 aero_state_total_num_conc = aero_state_total_num_conc & 00963 + aero_state_particle_num_conc(aero_state, & 00964 aero_state%apa%particle(i_part)) 00965 end do 00966 00967 end function aero_state_total_num_conc 00968 00969 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00970 00971 !> Does the same thing as aero_state_to_bin() but based on dry radius. 00972 subroutine aero_state_to_binned_dry(bin_grid, aero_data, aero_state, & 00973 aero_binned) 00974 00975 !> Bin grid. 00976 type(bin_grid_t), intent(in) :: bin_grid 00977 !> Aerosol data. 00978 type(aero_data_t), intent(in) :: aero_data 00979 !> Aerosol state. 00980 type(aero_state_t), intent(in) :: aero_state 00981 !> Binned distributions. 00982 type(aero_binned_t), intent(inout) :: aero_binned 00983 00984 integer :: i_part, i_bin 00985 type(aero_particle_t), pointer :: aero_particle 00986 00987 aero_binned%num_conc = 0d0 00988 aero_binned%vol_conc = 0d0 00989 do i_part = 1,aero_state%apa%n_part 00990 aero_particle => aero_state%apa%particle(i_part) 00991 i_bin = bin_grid_particle_in_bin(bin_grid, & 00992 aero_particle_solute_radius(aero_particle, aero_data)) 00993 if ((i_bin < 1) .or. (i_bin > bin_grid%n_bin)) then 00994 call warn_msg(503871022, "particle ID " & 00995 // trim(integer_to_string(aero_particle%id)) & 00996 // " outside of bin_grid, discarding") 00997 else 00998 aero_binned%vol_conc(i_bin,:) = aero_binned%vol_conc(i_bin,:) & 00999 + aero_particle%vol & 01000 * aero_weight_array_num_conc(aero_state%aero_weight, & 01001 aero_particle) / bin_grid%log_width 01002 aero_binned%num_conc(i_bin) = aero_binned%num_conc(i_bin) & 01003 + aero_weight_array_num_conc(aero_state%aero_weight, & 01004 aero_particle) / bin_grid%log_width 01005 end if 01006 end do 01007 01008 end subroutine aero_state_to_binned_dry 01009 01010 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01011 01012 !> Doubles number of particles in the given weight group. 01013 subroutine aero_state_double(aero_state, i_group) 01014 01015 !> Aerosol state. 01016 type(aero_state_t), intent(inout) :: aero_state 01017 !> Weight group to double. 01018 integer, intent(in) :: i_group 01019 01020 integer :: i_part 01021 type(aero_particle_t) :: aero_particle 01022 01023 call aero_particle_allocate(aero_particle) 01024 do i_part = 1,aero_state%apa%n_part 01025 if (aero_state%apa%particle(i_part)%weight_group == i_group) then 01026 call aero_particle_copy(aero_state%apa%particle(i_part), aero_particle) 01027 call aero_particle_new_id(aero_particle) 01028 call aero_state_add_particle(aero_state, aero_particle) 01029 end if 01030 end do 01031 call aero_particle_deallocate(aero_particle) 01032 aero_state%valid_sort = .false. 01033 call aero_weight_scale_comp_vol(aero_state%aero_weight(i_group), 2d0) 01034 01035 end subroutine aero_state_double 01036 01037 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01038 01039 !> Remove approximately half of the particles in the given weight group. 01040 subroutine aero_state_halve(aero_state, i_group) 01041 01042 !> Aerosol state. 01043 type(aero_state_t), intent(inout) :: aero_state 01044 !> Weight group to halve. 01045 integer, intent(in) :: i_group 01046 01047 integer :: i_part 01048 type(aero_info_t) :: aero_info 01049 01050 call aero_info_allocate(aero_info) 01051 do i_part = aero_state%apa%n_part,1,-1 01052 if (aero_state%apa%particle(i_part)%weight_group == i_group) then 01053 if (pmc_random() < 0.5d0) then 01054 aero_info%id = aero_state%apa%particle(i_part)%id 01055 aero_info%action = AERO_INFO_HALVED 01056 call aero_state_remove_particle_with_info(aero_state, i_part, & 01057 aero_info) 01058 end if 01059 end if 01060 end do 01061 call aero_info_deallocate(aero_info) 01062 call aero_weight_scale_comp_vol(aero_state%aero_weight(i_group), 0.5d0) 01063 01064 end subroutine aero_state_halve 01065 01066 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01067 01068 !> Double or halve the particle population in each weight group to 01069 !> maintain close to \c n_part_ideal particles per process, 01070 !> allocated equally amongst the weight groups. 01071 subroutine aero_state_rebalance(aero_state, allow_doubling, allow_halving, & 01072 initial_state_warning) 01073 01074 !> Aerosol state. 01075 type(aero_state_t), intent(inout) :: aero_state 01076 !> Whether to allow doubling of the population. 01077 logical, intent(in) :: allow_doubling 01078 !> Whether to allow halving of the population. 01079 logical, intent(in) :: allow_halving 01080 !> Whether to warn due to initial state doubling/halving. 01081 logical, intent(in) :: initial_state_warning 01082 01083 integer :: i_group, n_group, global_n_part 01084 01085 n_group = size(aero_state%aero_weight) 01086 01087 ! if we have less than half the maximum number of particles then 01088 ! double until we fill up the array 01089 if (allow_doubling) then 01090 do i_group = 1,n_group 01091 global_n_part & 01092 = aero_state_total_particles_all_procs(aero_state, i_group) 01093 do while ((real(global_n_part, kind=dp) 01094 < aero_state%n_part_ideal / real(n_group, kind=dp) / 2d0) 01095 .and. (global_n_part > 0)) 01096 if (initial_state_warning) then 01097 call warn_msg(716882783, "doubling particles in initial " & 01098 // "condition") 01099 end if 01100 call aero_state_double(aero_state, i_group) 01101 global_n_part & 01102 = aero_state_total_particles_all_procs(aero_state, i_group) 01103 end do 01104 end do 01105 end if 01106 01107 ! same for halving if we have too many particles 01108 if (allow_halving) then 01109 do i_group = 1,n_group 01110 do while (real(aero_state_total_particles_all_procs(aero_state, & i_group), kind=dp) 01111 > aero_state%n_part_ideal / real(n_group, kind=dp) * 2d0) 01112 if (initial_state_warning) then 01113 call warn_msg(661936373, & 01114 "halving particles in initial condition") 01115 end if 01116 call aero_state_halve(aero_state, i_group) 01117 end do 01118 end do 01119 end if 01120 01121 end subroutine aero_state_rebalance 01122 01123 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01124 01125 !> Scale the computational volume of the given group by the given 01126 !> ratio, altering particle number as necessary to preserve the 01127 !> number concentration. 01128 subroutine aero_state_scale_comp_vol(aero_state, i_group, comp_vol_ratio) 01129 01130 !> Aerosol state. 01131 type(aero_state_t), intent(inout) :: aero_state 01132 !> Group number. 01133 integer, intent(in) :: i_group 01134 !> Ratio of <tt>new_comp_vol / old_comp_vol</tt>. 01135 real(kind=dp), intent(in) :: comp_vol_ratio 01136 01137 real(kind=dp) :: ratio 01138 integer :: i_part, i_remove, n_remove, i_entry 01139 type(aero_info_t) :: aero_info 01140 01141 ! We could use the ratio > 1 case unconditionally, but that would 01142 ! have higher variance for the ratio < 1 case than the current 01143 ! scheme. 01144 01145 call aero_weight_scale_comp_vol(aero_state%aero_weight(i_group), & 01146 comp_vol_ratio) 01147 01148 if (aero_state%apa%n_part == 0) return 01149 01150 call aero_state_sort(aero_state) 01151 01152 if (comp_vol_ratio < 1d0) then 01153 n_remove = prob_round(comp_vol_ratio & 01154 * real(aero_state%aero_sorted%weight%inverse(i_group)%n_entry, kind=dp)) 01155 do i_remove = 1,n_remove 01156 i_entry = pmc_rand_int(aero_state%aero_sorted%weight%inverse(i_group)%n_entry) 01157 i_part = aero_state%aero_sorted%weight%inverse(i_group)%entry(i_entry) 01158 call aero_info_allocate(aero_info) 01159 aero_info%id = aero_state%apa%particle(i_part)%id 01160 aero_info%action = AERO_INFO_HALVED 01161 call aero_state_remove_particle(aero_state, i_part, .true., & 01162 aero_info) 01163 call aero_info_deallocate(aero_info) 01164 end do 01165 elseif (comp_vol_ratio > 1d0) then 01166 do i_entry = aero_state%aero_sorted%weight%inverse(i_group)%n_entry,1,-1 01167 i_part = aero_state%aero_sorted%weight%inverse(i_group)%entry(i_entry) 01168 call aero_state_dup_particle(aero_state, i_part, comp_vol_ratio) 01169 end do 01170 end if 01171 01172 end subroutine aero_state_scale_comp_vol 01173 01174 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01175 01176 !> Mix the aero_states between all processes. Currently uses a 01177 !> simple all-to-all diffusion. 01178 subroutine aero_state_mix(aero_state, del_t, mix_timescale, & 01179 aero_data, specify_prob_transfer) 01180 01181 !> Aerosol state. 01182 type(aero_state_t), intent(inout) :: aero_state 01183 !> Timestep (s). 01184 real(kind=dp), intent(in) :: del_t 01185 !> Mixing timescale (s). 01186 real(kind=dp), intent(in) :: mix_timescale 01187 !> Aero data values. 01188 type(aero_data_t), intent(in) :: aero_data 01189 !> Transfer probability of each particle (0 means no mixing, 1 01190 !> means total mixing). 01191 real(kind=dp), optional, intent(in) :: specify_prob_transfer 01192 01193 #ifdef PMC_USE_MPI 01194 integer :: rank, n_proc, i_proc, ierr 01195 integer :: buffer_size, buffer_size_check 01196 character, allocatable :: buffer(:) 01197 type(aero_state_t), allocatable :: aero_state_sends(:) 01198 type(aero_state_t), allocatable :: aero_state_recvs(:) 01199 real(kind=dp) :: prob_transfer, prob_not_transferred 01200 real(kind=dp) :: prob_transfer_given_not_transferred 01201 01202 ! process information 01203 rank = pmc_mpi_rank() 01204 n_proc = pmc_mpi_size() 01205 if (n_proc == 1) then 01206 ! buffer allocation below fails if n_proc == 1 01207 ! so bail out early (nothing to mix anyway) 01208 return 01209 end if 01210 01211 ! allocate aero_state arrays 01212 allocate(aero_state_sends(n_proc)) 01213 allocate(aero_state_recvs(n_proc)) 01214 do i_proc = 0,(n_proc - 1) 01215 call aero_state_allocate(aero_state_sends(i_proc + 1)) 01216 call aero_state_allocate(aero_state_recvs(i_proc + 1)) 01217 end do 01218 01219 ! compute the transfer probability 01220 if (present(specify_prob_transfer)) then 01221 prob_transfer = specify_prob_transfer / real(n_proc, kind=dp) 01222 else 01223 prob_transfer = (1d0 - exp(- del_t / mix_timescale)) & 01224 / real(n_proc, kind=dp) 01225 end if 01226 01227 ! extract particles to send 01228 prob_not_transferred = 1d0 01229 do i_proc = 0,(n_proc - 1) 01230 if (i_proc /= rank) then 01231 ! because we are doing sequential sampling from the aero_state 01232 ! we need to scale up the later transfer probabilities, because 01233 ! the later particles are being transferred conditioned on the 01234 ! fact that they did not transfer already 01235 prob_transfer_given_not_transferred = prob_transfer & 01236 / prob_not_transferred 01237 call aero_state_sample(aero_state, & 01238 aero_state_sends(i_proc + 1), & 01239 prob_transfer_given_not_transferred, AERO_INFO_NONE) 01240 prob_not_transferred = prob_not_transferred - prob_transfer 01241 end if 01242 end do 01243 01244 ! exchange the particles 01245 call aero_state_mpi_alltoall(aero_state_sends, aero_state_recvs) 01246 01247 ! process the received particles 01248 do i_proc = 0,(n_proc - 1) 01249 if (i_proc /= rank) then 01250 call aero_state_add(aero_state, aero_state_recvs(i_proc + 1)) 01251 end if 01252 end do 01253 01254 ! cleanup 01255 do i_proc = 0,(n_proc - 1) 01256 call aero_state_deallocate(aero_state_sends(i_proc + 1)) 01257 call aero_state_deallocate(aero_state_recvs(i_proc + 1)) 01258 end do 01259 deallocate(aero_state_sends) 01260 deallocate(aero_state_recvs) 01261 #endif 01262 01263 end subroutine aero_state_mix 01264 01265 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01266 01267 !> Do an MPI all-to-all transfer of aerosol states. 01268 !! 01269 !! States without particles are not sent. 01270 subroutine aero_state_mpi_alltoall(send, recv) 01271 01272 !> Array of aero_states to send (one per process). 01273 type(aero_state_t), intent(in) :: send(:) 01274 !> Array of aero_states to receives (one per process). 01275 type(aero_state_t), intent(inout) :: recv(size(send)) 01276 01277 #ifdef PMC_USE_MPI 01278 character, allocatable :: sendbuf(:), recvbuf(:) 01279 integer :: sendcounts(size(send)), sdispls(size(send)) 01280 integer :: recvcounts(size(send)), rdispls(size(send)) 01281 integer :: i_proc, position, old_position, max_sendbuf_size, ierr 01282 01283 call assert(978709842, size(send) == pmc_mpi_size()) 01284 01285 max_sendbuf_size = 0 01286 do i_proc = 1,pmc_mpi_size() 01287 if (aero_state_total_particles(send(i_proc)) > 0) then 01288 max_sendbuf_size = max_sendbuf_size & 01289 + pmc_mpi_pack_size_aero_state(send(i_proc)) 01290 end if 01291 end do 01292 01293 allocate(sendbuf(max_sendbuf_size)) 01294 01295 position = 0 01296 do i_proc = 1,pmc_mpi_size() 01297 old_position = position 01298 if (aero_state_total_particles(send(i_proc)) > 0) then 01299 call pmc_mpi_pack_aero_state(sendbuf, position, send(i_proc)) 01300 end if 01301 sendcounts(i_proc) = position - old_position 01302 end do 01303 call assert(393267406, position <= max_sendbuf_size) 01304 01305 call pmc_mpi_alltoall_integer(sendcounts, recvcounts) 01306 allocate(recvbuf(sum(recvcounts))) 01307 01308 sdispls(1) = 0 01309 rdispls(1) = 0 01310 do i_proc = 2,pmc_mpi_size() 01311 sdispls(i_proc) = sdispls(i_proc - 1) + sendcounts(i_proc - 1) 01312 rdispls(i_proc) = rdispls(i_proc - 1) + recvcounts(i_proc - 1) 01313 end do 01314 01315 call mpi_alltoallv(sendbuf, sendcounts, sdispls, MPI_CHARACTER, recvbuf, & 01316 recvcounts, rdispls, MPI_CHARACTER, MPI_COMM_WORLD, ierr) 01317 call pmc_mpi_check_ierr(ierr) 01318 01319 position = 0 01320 do i_proc = 1,pmc_mpi_size() 01321 call assert(189739257, position == rdispls(i_proc)) 01322 if (recvcounts(i_proc) > 0) then 01323 call pmc_mpi_unpack_aero_state(recvbuf, position, recv(i_proc)) 01324 end if 01325 end do 01326 01327 deallocate(sendbuf) 01328 deallocate(recvbuf) 01329 #endif 01330 01331 end subroutine aero_state_mpi_alltoall 01332 01333 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01334 01335 !> Set each aerosol particle to have its original total volume, but 01336 !> species volume ratios given by the total species volume ratio 01337 !> within each bin. This preserves the (weighted) total species 01338 !> volume per bin as well as per-particle total volumes. 01339 subroutine aero_state_bin_average_comp(aero_state, bin_grid, aero_data, & 01340 dry_volume) 01341 01342 !> Aerosol state to average. 01343 type(aero_state_t), intent(inout) :: aero_state 01344 !> Bin grid to average within. 01345 type(bin_grid_t), intent(in) :: bin_grid 01346 !> Aerosol data. 01347 type(aero_data_t), intent(in) :: aero_data 01348 !> Whether to use dry volume (rather than wet). 01349 logical, intent(in) :: dry_volume 01350 01351 real(kind=dp) :: species_volume_conc(aero_data%n_spec) 01352 real(kind=dp) :: total_volume_conc, particle_volume, num_conc 01353 integer :: i_bin, i_entry, i_part, i_spec 01354 type(aero_particle_t), pointer :: aero_particle 01355 01356 call aero_state_sort(aero_state, bin_grid) 01357 01358 do i_bin = 1,bin_grid%n_bin 01359 species_volume_conc = 0d0 01360 total_volume_conc = 0d0 01361 do i_entry = 1,aero_state%aero_sorted%size%inverse(i_bin)%n_entry 01362 i_part = aero_state%aero_sorted%size%inverse(i_bin)%entry(i_entry) 01363 aero_particle => aero_state%apa%particle(i_part) 01364 num_conc = aero_weight_array_num_conc(aero_state%aero_weight, & 01365 aero_particle) 01366 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01367 aero_data, dry_volume) 01368 species_volume_conc = species_volume_conc & 01369 + num_conc * aero_particle%vol 01370 total_volume_conc = total_volume_conc + num_conc * particle_volume 01371 end do 01372 do i_entry = 1,aero_state%aero_sorted%size%inverse(i_bin)%n_entry 01373 i_part = aero_state%aero_sorted%size%inverse(i_bin)%entry(i_entry) 01374 aero_particle => aero_state%apa%particle(i_part) 01375 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01376 aero_data, dry_volume) 01377 aero_particle%vol = particle_volume * species_volume_conc & 01378 / total_volume_conc 01379 if (dry_volume .and. (aero_data%i_water > 0)) then 01380 ! set water to zero if we are doing dry volume averaging 01381 aero_particle%vol(aero_data%i_water) = 0d0 01382 end if 01383 end do 01384 end do 01385 01386 end subroutine aero_state_bin_average_comp 01387 01388 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01389 01390 !> Set each aerosol particle to have its original species ratios, 01391 !> but total volume given by the average volume of all particles 01392 !> within each bin. 01393 !! 01394 !! This does not preserve the total species volume 01395 !! per bin. If the \c bin_center parameter is \c .true. then the 01396 !! particles in each bin are set to have the bin center volume, 01397 !! rather than the average volume of the particles in that bin. 01398 !! 01399 !! If the weighting function is not constant (AERO_WEIGHT_TYPE_NONE) 01400 !! then the averaging can be performed in either a number-preserving 01401 !! way or in a volume-preserving way. The volume-preserving way does 01402 !! not preserve species volume ratios in gernal, but will do so if 01403 !! the particle population has already been composition-averaged. 01404 subroutine aero_state_bin_average_size(aero_state, bin_grid, aero_data, & 01405 dry_volume, bin_center, preserve_number) 01406 01407 !> Aerosol state to average. 01408 type(aero_state_t), intent(inout) :: aero_state 01409 !> Bin grid to average within. 01410 type(bin_grid_t), intent(in) :: bin_grid 01411 !> Aerosol data. 01412 type(aero_data_t), intent(in) :: aero_data 01413 !> Whether to use dry volume (rather than wet). 01414 logical, intent(in) :: dry_volume 01415 !> Whether to assign the bin center volume (rather than the average 01416 !> volume). 01417 logical, intent(in) :: bin_center 01418 !> Whether to use the number-preserving scheme (otherwise will use 01419 !> the volume-preserving scheme). This parameter has no effect if 01420 !> \c bin_center is \c .true. 01421 logical, intent(in) :: preserve_number 01422 01423 real(kind=dp) :: total_volume_conc, particle_volume 01424 real(kind=dp) :: new_particle_volume, num_conc, total_num_conc 01425 real(kind=dp) :: lower_volume, upper_volume, center_volume 01426 real(kind=dp) :: lower_function, upper_function, center_function 01427 integer :: i_bin, i_entry, i_part, i_bisect, n_part 01428 logical :: monotone_increasing, monotone_decreasing 01429 type(aero_particle_t), pointer :: aero_particle 01430 01431 call aero_state_sort(aero_state, bin_grid) 01432 01433 do i_bin = 1,bin_grid%n_bin 01434 if (aero_state%aero_sorted%size%inverse(i_bin)%n_entry == 0) then 01435 cycle 01436 end if 01437 01438 n_part = aero_state%aero_sorted%size%inverse(i_bin)%n_entry 01439 total_num_conc = 0d0 01440 total_volume_conc = 0d0 01441 do i_entry = 1,aero_state%aero_sorted%size%inverse(i_bin)%n_entry 01442 i_part = aero_state%aero_sorted%size%inverse(i_bin)%entry(i_entry) 01443 aero_particle => aero_state%apa%particle(i_part) 01444 num_conc = aero_weight_array_num_conc(aero_state%aero_weight, & 01445 aero_particle) 01446 total_num_conc = total_num_conc + num_conc 01447 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01448 aero_data, dry_volume) 01449 total_volume_conc = total_volume_conc & 01450 + num_conc * particle_volume 01451 end do 01452 01453 ! determine the new_particle_volume for all particles in this bin 01454 if (bin_center) then 01455 new_particle_volume = rad2vol(bin_grid%center_radius(i_bin)) 01456 elseif (aero_weight_array_check_flat(aero_state%aero_weight)) then 01457 num_conc & ! any radius will have the same num_conc 01458 = aero_weight_array_num_conc_at_radius(aero_state%aero_weight, & 01459 1d0) 01460 new_particle_volume = total_volume_conc / num_conc & 01461 / real(aero_state%aero_sorted%size%inverse(i_bin)%n_entry, kind=dp) 01462 elseif (preserve_number) then 01463 ! number-preserving scheme: Solve the implicit equation: 01464 ! n_part * W(new_vol) = total_num_conc 01465 ! 01466 ! We assume that the weighting function is strictly monotone 01467 ! so this equation has a unique solution and the solution 01468 ! lies between the min and max particle volumes. We use 01469 ! bisection as this doesn't really need to be fast, just 01470 ! robust. 01471 01472 call aero_weight_array_check_monotonicity(aero_state%aero_weight, & 01473 monotone_increasing, monotone_decreasing) 01474 call assert_msg(214077200, & 01475 monotone_increasing .or. monotone_decreasing, & 01476 "monotone weight function required for averaging") 01477 01478 ! initialize to min and max particle volumes 01479 do i_entry = 1,n_part 01480 i_part = aero_state%aero_sorted%size%inverse(i_bin)%entry(i_entry) 01481 aero_particle => aero_state%apa%particle(i_part) 01482 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01483 aero_data, dry_volume) 01484 if (i_part == 1) then 01485 lower_volume = particle_volume 01486 upper_volume = particle_volume 01487 end if 01488 lower_volume = min(lower_volume, particle_volume) 01489 upper_volume = max(upper_volume, particle_volume) 01490 end do 01491 01492 lower_function = real(n_part, kind=dp) 01493 * aero_weight_array_num_conc_at_radius( 01494 aero_state%aero_weight, vol2rad(lower_volume)) - total_num_conc 01495 upper_function = real(n_part, kind=dp) 01496 * aero_weight_array_num_conc_at_radius( 01497 aero_state%aero_weight, vol2rad(upper_volume)) - total_num_conc 01498 01499 ! do 50 rounds of bisection (2^50 = 10^15) 01500 do i_bisect = 1,50 01501 center_volume = (lower_volume + upper_volume) / 2d0 01502 center_function = real(n_part, kind=dp) 01503 * aero_weight_array_num_conc_at_radius( 01504 aero_state%aero_weight, vol2rad(center_volume)) 01505 - total_num_conc 01506 if ((lower_function > 0d0 .and. center_function > 0d0) & 01507 .or. (lower_function < 0d0 .and. center_function < 0d0)) & 01508 then 01509 lower_volume = center_volume 01510 lower_function = center_function 01511 else 01512 upper_volume = center_volume 01513 upper_function = center_function 01514 end if 01515 end do 01516 01517 new_particle_volume = center_volume 01518 else 01519 ! volume-preserving scheme: Solve the implicit equation: 01520 ! n_part * W(new_vol) * new_vol = total_volume_conc 01521 ! 01522 ! We assume that the weighting function is strictly monotone 01523 ! so this equation has a unique solution and the solution 01524 ! lies between the min and max particle volumes. We use 01525 ! bisection as this doesn't really need to be fast, just 01526 ! robust. 01527 01528 call aero_weight_array_check_monotonicity(aero_state%aero_weight, & 01529 monotone_increasing, monotone_decreasing) 01530 call assert_msg(483078128, & 01531 monotone_increasing .or. monotone_decreasing, & 01532 "monotone weight function required for averaging") 01533 01534 ! initialize to min and max particle volumes 01535 do i_entry = 1,n_part 01536 i_part = aero_state%aero_sorted%size%inverse(i_bin)%entry(i_entry) 01537 aero_particle => aero_state%apa%particle(i_part) 01538 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01539 aero_data, dry_volume) 01540 if (i_part == 1) then 01541 lower_volume = particle_volume 01542 upper_volume = particle_volume 01543 end if 01544 lower_volume = min(lower_volume, particle_volume) 01545 upper_volume = max(upper_volume, particle_volume) 01546 end do 01547 01548 lower_function = real(n_part, kind=dp) 01549 * aero_weight_array_num_conc_at_radius( 01550 aero_state%aero_weight, vol2rad(lower_volume)) 01551 * lower_volume - total_volume_conc 01552 upper_function = real(n_part, kind=dp) 01553 * aero_weight_array_num_conc_at_radius( 01554 aero_state%aero_weight, vol2rad(upper_volume)) 01555 * upper_volume - total_volume_conc 01556 01557 ! do 50 rounds of bisection (2^50 = 10^15) 01558 do i_bisect = 1,50 01559 center_volume = (lower_volume + upper_volume) / 2d0 01560 center_function = real(n_part, kind=dp) 01561 * aero_weight_array_num_conc_at_radius( 01562 aero_state%aero_weight, vol2rad(center_volume)) 01563 * center_volume - total_volume_conc 01564 if ((lower_function > 0d0 .and. center_function > 0d0) & 01565 .or. (lower_function < 0d0 .and. center_function < 0d0)) & 01566 then 01567 lower_volume = center_volume 01568 lower_function = center_function 01569 else 01570 upper_volume = center_volume 01571 upper_function = center_function 01572 end if 01573 end do 01574 01575 new_particle_volume = center_volume 01576 end if 01577 01578 do i_entry = 1,n_part 01579 i_part = aero_state%aero_sorted%size%inverse(i_bin)%entry(i_entry) 01580 aero_particle => aero_state%apa%particle(i_part) 01581 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01582 aero_data, dry_volume) 01583 aero_particle%vol = aero_particle%vol / particle_volume & 01584 * new_particle_volume 01585 if (dry_volume .and. (aero_data%i_water > 0)) then 01586 ! set water to zero if we are doing dry volume averaging 01587 aero_particle%vol(aero_data%i_water) = 0d0 01588 end if 01589 end do 01590 end do 01591 01592 end subroutine aero_state_bin_average_size 01593 01594 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01595 01596 !> Make all particles dry (water set to zero). 01597 subroutine aero_state_make_dry(aero_state, aero_data) 01598 01599 !> Aerosol state to make dry. 01600 type(aero_state_t), intent(inout) :: aero_state 01601 !> Aerosol data. 01602 type(aero_data_t), intent(in) :: aero_data 01603 01604 integer :: i_part 01605 01606 if (aero_data%i_water > 0) then 01607 do i_part = 1,aero_state%apa%n_part 01608 aero_state%apa%particle(i_part)%vol(aero_data%i_water) = 0d0 01609 end do 01610 aero_state%valid_sort = .false. 01611 end if 01612 01613 end subroutine aero_state_make_dry 01614 01615 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01616 01617 !> Determines the number of bytes required to pack the given value. 01618 integer function pmc_mpi_pack_size_aero_state(val) 01619 01620 !> Value to pack. 01621 type(aero_state_t), intent(in) :: val 01622 01623 integer :: total_size, i_group 01624 01625 total_size = 0 01626 total_size = total_size + pmc_mpi_pack_size_apa(val%apa) 01627 total_size = total_size + pmc_mpi_pack_size_integer(size(val%aero_weight)) 01628 do i_group = 1,size(val%aero_weight) 01629 total_size = total_size & 01630 + pmc_mpi_pack_size_aero_weight(val%aero_weight(i_group)) 01631 end do 01632 total_size = total_size + pmc_mpi_pack_size_real(val%n_part_ideal) 01633 total_size = total_size + pmc_mpi_pack_size_aia(val%aero_info_array) 01634 pmc_mpi_pack_size_aero_state = total_size 01635 01636 end function pmc_mpi_pack_size_aero_state 01637 01638 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01639 01640 !> Packs the given value into the buffer, advancing position. 01641 subroutine pmc_mpi_pack_aero_state(buffer, position, val) 01642 01643 !> Memory buffer. 01644 character, intent(inout) :: buffer(:) 01645 !> Current buffer position. 01646 integer, intent(inout) :: position 01647 !> Value to pack. 01648 type(aero_state_t), intent(in) :: val 01649 01650 #ifdef PMC_USE_MPI 01651 integer :: prev_position, i_group 01652 01653 prev_position = position 01654 call pmc_mpi_pack_aero_particle_array(buffer, position, val%apa) 01655 call pmc_mpi_pack_integer(buffer, position, size(val%aero_weight)) 01656 do i_group = 1,size(val%aero_weight) 01657 call pmc_mpi_pack_aero_weight(buffer, position, & 01658 val%aero_weight(i_group)) 01659 end do 01660 call pmc_mpi_pack_real(buffer, position, val%n_part_ideal) 01661 call pmc_mpi_pack_aero_info_array(buffer, position, val%aero_info_array) 01662 call assert(850997402, & 01663 position - prev_position <= pmc_mpi_pack_size_aero_state(val)) 01664 #endif 01665 01666 end subroutine pmc_mpi_pack_aero_state 01667 01668 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01669 01670 !> Unpacks the given value from the buffer, advancing position. 01671 subroutine pmc_mpi_unpack_aero_state(buffer, position, val) 01672 01673 !> Memory buffer. 01674 character, intent(inout) :: buffer(:) 01675 !> Current buffer position. 01676 integer, intent(inout) :: position 01677 !> Value to pack. 01678 type(aero_state_t), intent(inout) :: val 01679 01680 #ifdef PMC_USE_MPI 01681 integer :: prev_position, i_group, n_group 01682 01683 val%valid_sort = .false. 01684 prev_position = position 01685 call pmc_mpi_unpack_aero_particle_array(buffer, position, val%apa) 01686 call pmc_mpi_unpack_integer(buffer, position, n_group) 01687 call aero_weight_deallocate(val%aero_weight) 01688 deallocate(val%aero_weight) 01689 allocate(val%aero_weight(n_group)) 01690 call aero_weight_allocate(val%aero_weight) 01691 do i_group = 1,n_group 01692 call pmc_mpi_unpack_aero_weight(buffer, position, & 01693 val%aero_weight(i_group)) 01694 end do 01695 call pmc_mpi_unpack_real(buffer, position, val%n_part_ideal) 01696 call pmc_mpi_unpack_aero_info_array(buffer, position, val%aero_info_array) 01697 call assert(132104747, & 01698 position - prev_position <= pmc_mpi_pack_size_aero_state(val)) 01699 #endif 01700 01701 end subroutine pmc_mpi_unpack_aero_state 01702 01703 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01704 01705 !> Gathers data from all processes into one aero_state on process 0. 01706 subroutine aero_state_mpi_gather(aero_state, aero_state_total) 01707 01708 !> Local aero_state. 01709 type(aero_state_t), intent(in) :: aero_state 01710 !> Centralized aero_state (only on process 0). 01711 type(aero_state_t), intent(inout) :: aero_state_total 01712 01713 #ifdef PMC_USE_MPI 01714 type(aero_state_t) :: aero_state_transfer 01715 integer :: n_proc, ierr, status(MPI_STATUS_SIZE) 01716 integer :: buffer_size, max_buffer_size, i_proc, position 01717 character, allocatable :: buffer(:) 01718 #endif 01719 01720 if (pmc_mpi_rank() == 0) then 01721 call aero_state_copy(aero_state, aero_state_total) 01722 end if 01723 01724 #ifdef PMC_USE_MPI 01725 01726 if (pmc_mpi_rank() /= 0) then 01727 ! send data from remote processes 01728 max_buffer_size = 0 01729 max_buffer_size = max_buffer_size & 01730 + pmc_mpi_pack_size_aero_state(aero_state) 01731 allocate(buffer(max_buffer_size)) 01732 position = 0 01733 call pmc_mpi_pack_aero_state(buffer, position, aero_state) 01734 call assert(542772170, position <= max_buffer_size) 01735 buffer_size = position 01736 call mpi_send(buffer, buffer_size, MPI_CHARACTER, 0, & 01737 AERO_STATE_TAG_GATHER, MPI_COMM_WORLD, ierr) 01738 call pmc_mpi_check_ierr(ierr) 01739 deallocate(buffer) 01740 else 01741 ! root process receives data from each remote proc 01742 n_proc = pmc_mpi_size() 01743 do i_proc = 1,(n_proc - 1) 01744 ! determine buffer size at root process 01745 call mpi_probe(i_proc, AERO_STATE_TAG_GATHER, MPI_COMM_WORLD, & 01746 status, ierr) 01747 call pmc_mpi_check_ierr(ierr) 01748 call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr) 01749 call pmc_mpi_check_ierr(ierr) 01750 01751 ! get buffer at root process 01752 allocate(buffer(buffer_size)) 01753 call mpi_recv(buffer, buffer_size, MPI_CHARACTER, i_proc, & 01754 AERO_STATE_TAG_GATHER, MPI_COMM_WORLD, status, ierr) 01755 01756 ! unpack it 01757 position = 0 01758 call aero_state_allocate(aero_state_transfer) 01759 call pmc_mpi_unpack_aero_state(buffer, position, & 01760 aero_state_transfer) 01761 call assert(518174881, position == buffer_size) 01762 deallocate(buffer) 01763 01764 call aero_state_add(aero_state_total, aero_state_transfer) 01765 01766 call aero_state_deallocate(aero_state_transfer) 01767 end do 01768 end if 01769 01770 #endif 01771 01772 end subroutine aero_state_mpi_gather 01773 01774 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01775 01776 !> Write the aero particle dimension to the given NetCDF file if it 01777 !> is not already present and in any case return the associated 01778 !> dimid. 01779 subroutine aero_state_netcdf_dim_aero_particle(aero_state, ncid, & 01780 dimid_aero_particle) 01781 01782 !> aero_state structure. 01783 type(aero_state_t), intent(in) :: aero_state 01784 !> NetCDF file ID, in data mode. 01785 integer, intent(in) :: ncid 01786 !> Dimid of the aero particle dimension. 01787 integer, intent(out) :: dimid_aero_particle 01788 01789 integer :: status, i_part 01790 integer :: varid_aero_particle 01791 integer :: aero_particle_centers(aero_state%apa%n_part) 01792 01793 ! try to get the dimension ID 01794 status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle) 01795 if (status == NF90_NOERR) return 01796 if (status /= NF90_EBADDIM) call pmc_nc_check(status) 01797 01798 ! dimension not defined, so define now define it 01799 call pmc_nc_check(nf90_redef(ncid)) 01800 01801 call pmc_nc_check(nf90_def_dim(ncid, "aero_particle", & 01802 aero_state%apa%n_part, dimid_aero_particle)) 01803 01804 call pmc_nc_check(nf90_enddef(ncid)) 01805 01806 do i_part = 1,aero_state%apa%n_part 01807 aero_particle_centers(i_part) = i_part 01808 end do 01809 call pmc_nc_write_integer_1d(ncid, aero_particle_centers, & 01810 "aero_particle", (/ dimid_aero_particle /), & 01811 description="dummy dimension variable (no useful value)") 01812 01813 end subroutine aero_state_netcdf_dim_aero_particle 01814 01815 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01816 01817 !> Write the aero removed dimension to the given NetCDF file if it 01818 !> is not already present and in any case return the associated 01819 !> dimid. 01820 subroutine aero_state_netcdf_dim_aero_removed(aero_state, ncid, & 01821 dimid_aero_removed) 01822 01823 !> aero_state structure. 01824 type(aero_state_t), intent(in) :: aero_state 01825 !> NetCDF file ID, in data mode. 01826 integer, intent(in) :: ncid 01827 !> Dimid of the aero removed dimension. 01828 integer, intent(out) :: dimid_aero_removed 01829 01830 integer :: status, i_remove, dim_size 01831 integer :: varid_aero_removed 01832 integer :: aero_removed_centers(max(aero_state%aero_info_array%n_item,1)) 01833 01834 ! try to get the dimension ID 01835 status = nf90_inq_dimid(ncid, "aero_removed", dimid_aero_removed) 01836 if (status == NF90_NOERR) return 01837 if (status /= NF90_EBADDIM) call pmc_nc_check(status) 01838 01839 ! dimension not defined, so define now define it 01840 call pmc_nc_check(nf90_redef(ncid)) 01841 01842 dim_size = max(aero_state%aero_info_array%n_item, 1) 01843 call pmc_nc_check(nf90_def_dim(ncid, "aero_removed", & 01844 dim_size, dimid_aero_removed)) 01845 01846 call pmc_nc_check(nf90_enddef(ncid)) 01847 01848 do i_remove = 1,dim_size 01849 aero_removed_centers(i_remove) = i_remove 01850 end do 01851 call pmc_nc_write_integer_1d(ncid, aero_removed_centers, & 01852 "aero_removed", (/ dimid_aero_removed /), & 01853 description="dummy dimension variable (no useful value)") 01854 01855 end subroutine aero_state_netcdf_dim_aero_removed 01856 01857 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01858 01859 !> Write full state. 01860 subroutine aero_state_output_netcdf(aero_state, ncid, aero_data, & 01861 record_removals, record_optical) 01862 01863 !> aero_state to write. 01864 type(aero_state_t), intent(in) :: aero_state 01865 !> NetCDF file ID, in data mode. 01866 integer, intent(in) :: ncid 01867 !> aero_data structure. 01868 type(aero_data_t), intent(in) :: aero_data 01869 !> Whether to output particle removal info. 01870 logical, intent(in) :: record_removals 01871 !> Whether to output aerosol optical properties. 01872 logical, intent(in) :: record_optical 01873 01874 integer :: dimid_aero_particle, dimid_aero_species, dimid_aero_source 01875 integer :: dimid_aero_removed 01876 integer :: i_part, i_remove 01877 type(aero_particle_t), pointer :: particle 01878 real(kind=dp) :: aero_particle_mass(aero_state%apa%n_part, aero_data%n_spec) 01879 integer :: aero_n_orig_part(aero_state%apa%n_part, aero_data%n_source) 01880 integer :: aero_weight_group(aero_state%apa%n_part) 01881 real(kind=dp) :: aero_absorb_cross_sect(aero_state%apa%n_part) 01882 real(kind=dp) :: aero_scatter_cross_sect(aero_state%apa%n_part) 01883 real(kind=dp) :: aero_asymmetry(aero_state%apa%n_part) 01884 real(kind=dp) :: aero_refract_shell_real(aero_state%apa%n_part) 01885 real(kind=dp) :: aero_refract_shell_imag(aero_state%apa%n_part) 01886 real(kind=dp) :: aero_refract_core_real(aero_state%apa%n_part) 01887 real(kind=dp) :: aero_refract_core_imag(aero_state%apa%n_part) 01888 real(kind=dp) :: aero_core_vol(aero_state%apa%n_part) 01889 integer :: aero_water_hyst_leg(aero_state%apa%n_part) 01890 real(kind=dp) :: aero_comp_vol(aero_state%apa%n_part) 01891 integer :: aero_id(aero_state%apa%n_part) 01892 real(kind=dp) :: aero_least_create_time(aero_state%apa%n_part) 01893 real(kind=dp) :: aero_greatest_create_time(aero_state%apa%n_part) 01894 integer :: aero_removed_id(max(aero_state%aero_info_array%n_item,1)) 01895 integer :: aero_removed_action(max(aero_state%aero_info_array%n_item,1)) 01896 integer :: aero_removed_other_id(max(aero_state%aero_info_array%n_item,1)) 01897 01898 !> \page output_format_aero_state Output File Format: Aerosol Particle State 01899 !! 01900 !! The aerosol state consists of a set of individual aerosol 01901 !! particles, each with its own individual properties. The 01902 !! properties of all particles are stored in arrays, one per 01903 !! property. For example, <tt>aero_absorb_cross_sect(i)</tt> gives 01904 !! the absorption cross section of particle number \c i, while 01905 !! <tt>aero_particle_mass(i,s)</tt> gives the mass of species \c s 01906 !! in particle \c i. The aerosol species are described in \ref 01907 !! output_format_aero_data. 01908 !! 01909 !! Each aerosol particle \c i represents a volume of space known 01910 !! as the computational volume and given by 01911 !! <tt>aero_comp_vol(i)</tt>. Dividing a per-particle quantity by 01912 !! the respective computational volume gives the concentration of 01913 !! that quantity contributed by the particle. For example, summing 01914 !! <tt>aero_particle_mass(i,s) / aero_comp_vol(i)</tt> over all \c 01915 !! i gives the total concentration of species \c s in 01916 !! kg/m^3. Similarly, summing <tt>aero_absorb_cross_sect(i) / 01917 !! aero_comp_vol(i)</tt> over all \c i will give the concentration 01918 !! of scattering cross section in m^2/m^3. 01919 !! 01920 !! FIXME: the aero_weight is also output 01921 !! 01922 !! The aerosol state uses the \c aero_species NetCDF dimension as 01923 !! specified in the \ref output_format_aero_data section, as well 01924 !! as the NetCDF dimension: 01925 !! - \b aero_particle: number of aerosol particles 01926 !! 01927 !! The aerosol state NetCDF variables are: 01928 !! - \b aero_particle (dim \c aero_particle): dummy dimension variable 01929 !! (no useful value) 01930 !! - \b aero_particle_mass (unit kg, 01931 !! dim <tt>aero_particle x aero_species</tt>): constituent masses of 01932 !! each aerosol particle - <tt>aero_particle_mass(i,s)</tt> gives the 01933 !! mass of species \c s in particle \c i 01934 !! - \b aero_n_orig_part (dim <tt>aero_particle x 01935 !! aero_source</tt>): number of original particles from each 01936 !! source that formed each aerosol particle - 01937 !! <tt>aero_n_orig_part(i,s)</tt> is the number of particles 01938 !! from source \c s that contributed to particle \c i - when 01939 !! particle \c i first enters the simulation (by emissions, 01940 !! dilution, etc.) it has <tt>aero_n_orig_part(i,s) = 1</tt> 01941 !! for the source number \c s it came from (otherwise zero) 01942 !! and when two particles coagulate, their values of \c 01943 !! aero_n_orig_part are added (the number of coagulation 01944 !! events that formed each particle is thus 01945 !! <tt>sum(aero_n_orig_part(i,:)) - 1</tt>) 01946 !! - \b aero_weight_group (dim <tt>aero_particle</tt>): weight 01947 !! group number (1 to <tt>aero_weight</tt>) of each aerosol 01948 !! particle 01949 !! - \b aero_absorb_cross_sect (unit m^2, dim \c aero_particle): 01950 !! optical absorption cross sections of each aerosol particle 01951 !! - \b aero_scatter_cross_sect (unit m^2, dim \c aero_particle): 01952 !! optical scattering cross sections of each aerosol particle 01953 !! - \b aero_asymmetry (dimensionless, dim \c aero_particle): optical 01954 !! asymmetry parameters of each aerosol particle 01955 !! - \b aero_refract_shell_real (dimensionless, dim \c aero_particle): 01956 !! real part of the refractive indices of the shell of each 01957 !! aerosol particle 01958 !! - \b aero_refract_shell_imag (dimensionless, dim \c aero_particle): 01959 !! imaginary part of the refractive indices of the shell of each 01960 !! aerosol particle 01961 !! - \b aero_refract_core_real (dimensionless, dim \c aero_particle): 01962 !! real part of the refractive indices of the core of each 01963 !! aerosol particle 01964 !! - \b aero_refract_core_imag (dimensionless, dim \c aero_particle): 01965 !! imaginary part of the refractive indices of the core of each 01966 !! aerosol particle 01967 !! - \b aero_core_vol (unit m^3, dim \c aero_particle): volume of the 01968 !! optical cores of each aerosol particle 01969 !! - \b aero_water_hyst_leg (dim \c aero_particle): integers 01970 !! specifying which leg of the water hysteresis curve each 01971 !! particle is on, using the MOSAIC numbering convention 01972 !! - \b aero_comp_vol (unit m^3, dim \c aero_particle): computational 01973 !! volume containing each particle 01974 !! - \b aero_id (dim \c aero_particle): unique ID number of each 01975 !! aerosol particle 01976 !! - \b aero_least_create_time (unit s, dim \c aero_particle): least 01977 !! (earliest) creation time of any original constituent particles 01978 !! that coagulated to form each particle, measured from the start 01979 !! of the simulation - a particle is said to be created when it 01980 !! first enters the simulation (by emissions, dilution, etc.) 01981 !! - \b aero_greatest_create_time (unit s, dim \c 01982 !! aero_particle): greatest (latest) creation time of any 01983 !! original constituent particles that coagulated to form each 01984 !! particle, measured from the start of the simulation - a 01985 !! particle is said to be created when it first enters the 01986 !! simulation (by emissions, dilution, etc.) 01987 01988 call aero_weight_array_output_netcdf(aero_state%aero_weight, ncid) 01989 01990 call aero_data_netcdf_dim_aero_species(aero_data, ncid, & 01991 dimid_aero_species) 01992 call aero_data_netcdf_dim_aero_source(aero_data, ncid, & 01993 dimid_aero_source) 01994 01995 if (aero_state%apa%n_part > 0) then 01996 call aero_state_netcdf_dim_aero_particle(aero_state, ncid, & 01997 dimid_aero_particle) 01998 01999 ! FIXME: replace this loop with statements like 02000 ! aero_n_orig_part = aero_state%apa%particle%n_orig_part 02001 do i_part = 1,aero_state%apa%n_part 02002 particle => aero_state%apa%particle(i_part) 02003 aero_particle_mass(i_part, :) = particle%vol * aero_data%density 02004 aero_n_orig_part(i_part, :) = particle%n_orig_part 02005 aero_weight_group(i_part) = particle%weight_group 02006 aero_water_hyst_leg(i_part) = particle%water_hyst_leg 02007 aero_comp_vol(i_part) & 02008 = 1d0 / aero_state_particle_num_conc(aero_state, particle) 02009 aero_id(i_part) = particle%id 02010 aero_least_create_time(i_part) = particle%least_create_time 02011 aero_greatest_create_time(i_part) = particle%greatest_create_time 02012 if (record_optical) then 02013 aero_absorb_cross_sect(i_part) = particle%absorb_cross_sect 02014 aero_scatter_cross_sect(i_part) = particle%scatter_cross_sect 02015 aero_asymmetry(i_part) = particle%asymmetry 02016 aero_refract_shell_real(i_part) = real(particle%refract_shell) 02017 aero_refract_shell_imag(i_part) = aimag(particle%refract_shell) 02018 aero_refract_core_real(i_part) = real(particle%refract_core) 02019 aero_refract_core_imag(i_part) = aimag(particle%refract_core) 02020 aero_core_vol(i_part) = particle%core_vol 02021 end if 02022 end do 02023 call pmc_nc_write_real_2d(ncid, aero_particle_mass, & 02024 "aero_particle_mass", (/ dimid_aero_particle, & 02025 dimid_aero_species /), unit="kg", & 02026 long_name="constituent masses of each aerosol particle") 02027 call pmc_nc_write_integer_2d(ncid, aero_n_orig_part, & 02028 "aero_n_orig_part", (/ dimid_aero_particle, & 02029 dimid_aero_source /), & 02030 long_name="number of original constituent particles from " & 02031 // "each source that coagulated to form each aerosol particle") 02032 call pmc_nc_write_integer_1d(ncid, aero_weight_group, & 02033 "aero_weight_group", (/ dimid_aero_particle /), & 02034 long_name="weight group number of each aerosol particle") 02035 call pmc_nc_write_integer_1d(ncid, aero_water_hyst_leg, & 02036 "aero_water_hyst_leg", (/ dimid_aero_particle /), & 02037 long_name="leg of the water hysteresis curve leg of each "& 02038 // "aerosol particle") 02039 call pmc_nc_write_real_1d(ncid, aero_comp_vol, & 02040 "aero_comp_vol", (/ dimid_aero_particle /), unit="m^3", & 02041 long_name="computational volume containing each particle") 02042 call pmc_nc_write_integer_1d(ncid, aero_id, & 02043 "aero_id", (/ dimid_aero_particle /), & 02044 long_name="unique ID number of each aerosol particle") 02045 call pmc_nc_write_real_1d(ncid, aero_least_create_time, & 02046 "aero_least_create_time", (/ dimid_aero_particle /), unit="s", & 02047 long_name="least creation time of each aerosol particle", & 02048 description="least (earliest) creation time of any original " & 02049 // "constituent particles that coagulated to form each " & 02050 // "particle, measured from the start of the simulation") 02051 call pmc_nc_write_real_1d(ncid, aero_greatest_create_time, & 02052 "aero_greatest_create_time", (/ dimid_aero_particle /), & 02053 unit="s", & 02054 long_name="greatest creation time of each aerosol particle", & 02055 description="greatest (latest) creation time of any original " & 02056 // "constituent particles that coagulated to form each " & 02057 // "particle, measured from the start of the simulation") 02058 if (record_optical) then 02059 call pmc_nc_write_real_1d(ncid, aero_absorb_cross_sect, & 02060 "aero_absorb_cross_sect", (/ dimid_aero_particle /), & 02061 unit="m^2", & 02062 long_name="optical absorption cross sections of each " & 02063 // "aerosol particle") 02064 call pmc_nc_write_real_1d(ncid, aero_scatter_cross_sect, & 02065 "aero_scatter_cross_sect", (/ dimid_aero_particle /), & 02066 unit="m^2", & 02067 long_name="optical scattering cross sections of each " & 02068 // "aerosol particle") 02069 call pmc_nc_write_real_1d(ncid, aero_asymmetry, & 02070 "aero_asymmetry", (/ dimid_aero_particle /), unit="1", & 02071 long_name="optical asymmetry parameters of each " & 02072 // "aerosol particle") 02073 call pmc_nc_write_real_1d(ncid, aero_refract_shell_real, & 02074 "aero_refract_shell_real", (/ dimid_aero_particle /), & 02075 unit="1", & 02076 long_name="real part of the refractive indices of the " & 02077 // "shell of each aerosol particle") 02078 call pmc_nc_write_real_1d(ncid, aero_refract_shell_imag, & 02079 "aero_refract_shell_imag", (/ dimid_aero_particle /), & 02080 unit="1", & 02081 long_name="imaginary part of the refractive indices of " & 02082 // "the shell of each aerosol particle") 02083 call pmc_nc_write_real_1d(ncid, aero_refract_core_real, & 02084 "aero_refract_core_real", (/ dimid_aero_particle /), & 02085 unit="1", & 02086 long_name="real part of the refractive indices of the core " & 02087 // "of each aerosol particle") 02088 call pmc_nc_write_real_1d(ncid, aero_refract_core_imag, & 02089 "aero_refract_core_imag", (/ dimid_aero_particle /), & 02090 unit="1", & 02091 long_name="imaginary part of the refractive indices of " & 02092 // "the core of each aerosol particle") 02093 call pmc_nc_write_real_1d(ncid, aero_core_vol, & 02094 "aero_core_vol", (/ dimid_aero_particle /), unit="m^3", & 02095 long_name="volume of the optical cores of each " & 02096 // "aerosol particle") 02097 end if 02098 end if 02099 02100 !> \page output_format_aero_removed Output File Format: Aerosol Particle Removal Information 02101 !! 02102 !! When an aerosol particle is introduced into the simulation it 02103 !! is assigned a unique ID number. This ID number will persist 02104 !! over time, allowing tracking of a paticular particle's 02105 !! evolution. If the \c record_removals variable in the input spec 02106 !! file is \c yes, then the every time a particle is removed from 02107 !! the simulation its removal will be recorded in the removal 02108 !! information. 02109 !! 02110 !! The removal information written at timestep \c n contains 02111 !! information about every particle ID that is present at time 02112 !! <tt>(n - 1)</tt> but not present at time \c n. 02113 !! 02114 !! The removal information is always written in the output files, 02115 !! even if no particles were removed in the previous 02116 !! timestep. Unfortunately, NetCDF files cannot contain arrays of 02117 !! length 0. In the case of no particles being removed, the \c 02118 !! aero_removed dimension will be set to 1 and 02119 !! <tt>aero_removed_action(1)</tt> will be 0 (\c AERO_INFO_NONE). 02120 !! 02121 !! When two particles coagulate, the ID number of the combined 02122 !! particle will be the ID particle of the largest constituent, if 02123 !! possible (weighting functions can make this impossible to 02124 !! achieve). A given particle ID may thus be lost due to 02125 !! coagulation (if the resulting combined particle has a different 02126 !! ID), or the ID may be preserved (as the ID of the combined 02127 !! particle). Only if the ID is lost will the particle be recorded 02128 !! in the removal information, and in this case 02129 !! <tt>aero_removed_action(i)</tt> will be 2 (\c AERO_INFO_COAG) 02130 !! and <tt>aero_removed_other_id(i)</tt> will be the ID number of 02131 !! the combined particle. 02132 !! 02133 !! The aerosol removal information NetCDF dimensions are: 02134 !! - \b aero_removed: number of aerosol particles removed from the 02135 !! simulation during the previous timestep (or 1, as described 02136 !! above) 02137 !! 02138 !! The aerosol removal information NetCDF variables are: 02139 !! - \b aero_removed (dim \c aero_removed): dummy dimension variable 02140 !! (no useful value) 02141 !! - \b aero_removed_id (dim \c aero_removed): the ID number of each 02142 !! removed particle 02143 !! - \b aero_removed_action (dim \c aero_removed): the reasons for 02144 !! removal for each particle, with values: 02145 !! - 0 (\c AERO_INFO_NONE): no information (invalid entry) 02146 !! - 1 (\c AERO_INFO_DILUTION): particle was removed due to dilution 02147 !! with outside air 02148 !! - 2 (\c AERO_INFO_COAG): particle was removed due to coagulation 02149 !! - 3 (\c AERO_INFO_HALVED): particle was removed due to halving of 02150 !! the aerosol population 02151 !! - 4 (\c AERO_INFO_WEIGHT): particle was removed due to adjustments 02152 !! in the particle's weighting function 02153 !! - \b aero_removed_other_id (dim \c aero_removed): the ID number of 02154 !! the combined particle formed by coagulation, if the removal reason 02155 !! was coagulation (2, \c AERO_INFO_COAG). May be 0, if the new 02156 !! coagulated particle was not created due to weighting. 02157 02158 ! FIXME: move this to aero_info_array.F90, together with 02159 ! aero_state_netcdf_dim_aero_removed() ? 02160 if (record_removals) then 02161 call aero_state_netcdf_dim_aero_removed(aero_state, ncid, & 02162 dimid_aero_removed) 02163 if (aero_state%aero_info_array%n_item >= 1) then 02164 do i_remove = 1,aero_state%aero_info_array%n_item 02165 aero_removed_id(i_remove) = & 02166 aero_state%aero_info_array%aero_info(i_remove)%id 02167 aero_removed_action(i_remove) = & 02168 aero_state%aero_info_array%aero_info(i_remove)%action 02169 aero_removed_other_id(i_remove) = & 02170 aero_state%aero_info_array%aero_info(i_remove)%other_id 02171 end do 02172 else 02173 aero_removed_id(1) = 0 02174 aero_removed_action(1) = AERO_INFO_NONE 02175 aero_removed_other_id(1) = 0 02176 end if 02177 call pmc_nc_write_integer_1d(ncid, aero_removed_id, & 02178 "aero_removed_id", (/ dimid_aero_removed /), & 02179 long_name="ID of removed particles") 02180 call pmc_nc_write_integer_1d(ncid, aero_removed_action, & 02181 "aero_removed_action", (/ dimid_aero_removed /), & 02182 long_name="reason for particle removal", & 02183 description="valid is 0 (invalid entry), 1 (removed due to " & 02184 // "dilution), 2 (removed due to coagulation -- combined " & 02185 // "particle ID is in \c aero_removed_other_id), 3 (removed " & 02186 // "due to populating halving), or 4 (removed due to " & 02187 // "weighting changes") 02188 call pmc_nc_write_integer_1d(ncid, aero_removed_other_id, & 02189 "aero_removed_other_id", (/ dimid_aero_removed /), & 02190 long_name="ID of other particle involved in removal", & 02191 description="if <tt>aero_removed_action(i)</tt> is 2 " & 02192 // "(due to coagulation), then " & 02193 // "<tt>aero_removed_other_id(i)</tt> is the ID of the " & 02194 // "resulting combined particle, or 0 if the new particle " & 02195 // "was not created") 02196 end if 02197 02198 end subroutine aero_state_output_netcdf 02199 02200 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02201 02202 !> Read full state. 02203 subroutine aero_state_input_netcdf(aero_state, ncid, aero_data) 02204 02205 !> aero_state to read. 02206 type(aero_state_t), intent(inout) :: aero_state 02207 !> NetCDF file ID, in data mode. 02208 integer, intent(in) :: ncid 02209 !> aero_data structure. 02210 type(aero_data_t), intent(in) :: aero_data 02211 02212 integer :: dimid_aero_particle, dimid_aero_removed, n_info_item, n_part 02213 integer :: i_bin, i_part_in_bin, i_part, i_remove, status 02214 type(aero_particle_t) :: aero_particle 02215 character(len=1000) :: name 02216 02217 real(kind=dp), allocatable :: aero_particle_mass(:,:) 02218 integer, allocatable :: aero_n_orig_part(:,:) 02219 integer, allocatable :: aero_weight_group(:) 02220 real(kind=dp), allocatable :: aero_absorb_cross_sect(:) 02221 real(kind=dp), allocatable :: aero_scatter_cross_sect(:) 02222 real(kind=dp), allocatable :: aero_asymmetry(:) 02223 real(kind=dp), allocatable :: aero_refract_shell_real(:) 02224 real(kind=dp), allocatable :: aero_refract_shell_imag(:) 02225 real(kind=dp), allocatable :: aero_refract_core_real(:) 02226 real(kind=dp), allocatable :: aero_refract_core_imag(:) 02227 real(kind=dp), allocatable :: aero_core_vol(:) 02228 integer, allocatable :: aero_water_hyst_leg(:) 02229 real(kind=dp), allocatable :: aero_comp_vol(:) 02230 integer, allocatable :: aero_id(:) 02231 real(kind=dp), allocatable :: aero_least_create_time(:) 02232 real(kind=dp), allocatable :: aero_greatest_create_time(:) 02233 integer, allocatable :: aero_removed_id(:) 02234 integer, allocatable :: aero_removed_action(:) 02235 integer, allocatable :: aero_removed_other_id(:) 02236 02237 status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle) 02238 if (status == NF90_EBADDIM) then 02239 ! no aero_particle dimension means no particles present 02240 call aero_state_deallocate(aero_state) 02241 call aero_state_allocate(aero_state) 02242 call aero_weight_array_input_netcdf(aero_state%aero_weight, ncid) 02243 return 02244 end if 02245 call pmc_nc_check(status) 02246 call pmc_nc_check(nf90_Inquire_Dimension(ncid, dimid_aero_particle, & 02247 name, n_part)) 02248 02249 allocate(aero_particle_mass(n_part, aero_data%n_spec)) 02250 allocate(aero_n_orig_part(n_part, aero_data%n_source)) 02251 allocate(aero_weight_group(n_part)) 02252 allocate(aero_absorb_cross_sect(n_part)) 02253 allocate(aero_scatter_cross_sect(n_part)) 02254 allocate(aero_asymmetry(n_part)) 02255 allocate(aero_refract_shell_real(n_part)) 02256 allocate(aero_refract_shell_imag(n_part)) 02257 allocate(aero_refract_core_real(n_part)) 02258 allocate(aero_refract_core_imag(n_part)) 02259 allocate(aero_core_vol(n_part)) 02260 allocate(aero_water_hyst_leg(n_part)) 02261 allocate(aero_comp_vol(n_part)) 02262 allocate(aero_id(n_part)) 02263 allocate(aero_least_create_time(n_part)) 02264 allocate(aero_greatest_create_time(n_part)) 02265 02266 call pmc_nc_read_real_2d(ncid, aero_particle_mass, & 02267 "aero_particle_mass") 02268 call pmc_nc_read_integer_2d(ncid, aero_n_orig_part, & 02269 "aero_n_orig_part") 02270 call pmc_nc_read_integer_1d(ncid, aero_weight_group, & 02271 "aero_weight_group") 02272 call pmc_nc_read_real_1d(ncid, aero_absorb_cross_sect, & 02273 "aero_absorb_cross_sect", must_be_present=.false.) 02274 call pmc_nc_read_real_1d(ncid, aero_scatter_cross_sect, & 02275 "aero_scatter_cross_sect", must_be_present=.false.) 02276 call pmc_nc_read_real_1d(ncid, aero_asymmetry, & 02277 "aero_asymmetry", must_be_present=.false.) 02278 call pmc_nc_read_real_1d(ncid, aero_refract_shell_real, & 02279 "aero_refract_shell_real", must_be_present=.false.) 02280 call pmc_nc_read_real_1d(ncid, aero_refract_shell_imag, & 02281 "aero_refract_shell_imag", must_be_present=.false.) 02282 call pmc_nc_read_real_1d(ncid, aero_refract_core_real, & 02283 "aero_refract_core_real", must_be_present=.false.) 02284 call pmc_nc_read_real_1d(ncid, aero_refract_core_imag, & 02285 "aero_refract_core_imag", must_be_present=.false.) 02286 call pmc_nc_read_real_1d(ncid, aero_core_vol, & 02287 "aero_core_vol", must_be_present=.false.) 02288 call pmc_nc_read_integer_1d(ncid, aero_water_hyst_leg, & 02289 "aero_water_hyst_leg") 02290 call pmc_nc_read_real_1d(ncid, aero_comp_vol, & 02291 "aero_comp_vol") 02292 call pmc_nc_read_integer_1d(ncid, aero_id, & 02293 "aero_id") 02294 call pmc_nc_read_real_1d(ncid, aero_least_create_time, & 02295 "aero_least_create_time") 02296 call pmc_nc_read_real_1d(ncid, aero_greatest_create_time, & 02297 "aero_greatest_create_time") 02298 02299 call aero_state_deallocate(aero_state) 02300 call aero_state_allocate(aero_state) 02301 02302 call aero_weight_array_input_netcdf(aero_state%aero_weight, ncid) 02303 aero_state%n_part_ideal = 0d0 02304 02305 call aero_particle_allocate_size(aero_particle, aero_data%n_spec, & 02306 aero_data%n_source) 02307 do i_part = 1,n_part 02308 aero_particle%vol = aero_particle_mass(i_part, :) / aero_data%density 02309 aero_particle%n_orig_part = aero_n_orig_part(i_part, :) 02310 aero_particle%weight_group = aero_weight_group(i_part) 02311 aero_particle%absorb_cross_sect = aero_absorb_cross_sect(i_part) 02312 aero_particle%scatter_cross_sect = aero_scatter_cross_sect(i_part) 02313 aero_particle%asymmetry = aero_asymmetry(i_part) 02314 aero_particle%refract_shell = & 02315 cmplx(aero_refract_shell_real(i_part), & 02316 aero_refract_shell_imag(i_part), kind=dc) 02317 aero_particle%refract_core = cmplx(aero_refract_core_real(i_part), & 02318 aero_refract_core_imag(i_part), kind=dc) 02319 aero_particle%core_vol = aero_core_vol(i_part) 02320 aero_particle%water_hyst_leg = aero_water_hyst_leg(i_part) 02321 aero_particle%id = aero_id(i_part) 02322 aero_particle%least_create_time = aero_least_create_time(i_part) 02323 aero_particle%greatest_create_time = aero_greatest_create_time(i_part) 02324 02325 call assert(314368871, almost_equal(aero_comp_vol(i_part), & 02326 1d0 / aero_weight_array_num_conc(aero_state%aero_weight, & 02327 aero_particle))) 02328 02329 call aero_state_add_particle(aero_state, aero_particle) 02330 end do 02331 call aero_particle_deallocate(aero_particle) 02332 02333 deallocate(aero_particle_mass) 02334 deallocate(aero_n_orig_part) 02335 deallocate(aero_weight_group) 02336 deallocate(aero_absorb_cross_sect) 02337 deallocate(aero_scatter_cross_sect) 02338 deallocate(aero_asymmetry) 02339 deallocate(aero_refract_shell_real) 02340 deallocate(aero_refract_shell_imag) 02341 deallocate(aero_refract_core_real) 02342 deallocate(aero_refract_core_imag) 02343 deallocate(aero_core_vol) 02344 deallocate(aero_water_hyst_leg) 02345 deallocate(aero_comp_vol) 02346 deallocate(aero_id) 02347 deallocate(aero_least_create_time) 02348 deallocate(aero_greatest_create_time) 02349 02350 status = nf90_inq_dimid(ncid, "aero_removed", dimid_aero_removed) 02351 if ((status /= NF90_NOERR) .and. (status /= NF90_EBADDIM)) then 02352 call pmc_nc_check(status) 02353 end if 02354 if (status == NF90_NOERR) then 02355 call pmc_nc_check(nf90_Inquire_Dimension(ncid, dimid_aero_removed, & 02356 name, n_info_item)) 02357 02358 allocate(aero_removed_id(max(n_info_item,1))) 02359 allocate(aero_removed_action(max(n_info_item,1))) 02360 allocate(aero_removed_other_id(max(n_info_item,1))) 02361 02362 call pmc_nc_read_integer_1d(ncid, aero_removed_id, & 02363 "aero_removed_id") 02364 call pmc_nc_read_integer_1d(ncid, aero_removed_action, & 02365 "aero_removed_action") 02366 call pmc_nc_read_integer_1d(ncid, aero_removed_other_id, & 02367 "aero_removed_other_id") 02368 02369 if ((n_info_item > 1) .or. (aero_removed_id(1) /= 0)) then 02370 call aero_info_array_enlarge_to(aero_state%aero_info_array, & 02371 n_info_item) 02372 do i_remove = 1,n_info_item 02373 aero_state%aero_info_array%aero_info(i_remove)%id & 02374 = aero_removed_id(i_remove) 02375 aero_state%aero_info_array%aero_info(i_remove)%action & 02376 = aero_removed_action(i_remove) 02377 aero_state%aero_info_array%aero_info(i_remove)%other_id & 02378 = aero_removed_other_id(i_remove) 02379 end do 02380 end if 02381 02382 deallocate(aero_removed_id) 02383 deallocate(aero_removed_action) 02384 deallocate(aero_removed_other_id) 02385 end if 02386 02387 end subroutine aero_state_input_netcdf 02388 02389 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02390 02391 !> Sorts the particles if necessary. 02392 subroutine aero_state_sort(aero_state, bin_grid, all_procs_same) 02393 02394 !> Aerosol state to sort. 02395 type(aero_state_t), intent(inout) :: aero_state 02396 !> Bin grid. 02397 type(bin_grid_t), optional, intent(in) :: bin_grid 02398 !> Whether all processors should use the same bin grid. 02399 logical, optional, intent(in) :: all_procs_same 02400 02401 call aero_sorted_remake_if_needed(aero_state%aero_sorted, aero_state%apa, & 02402 aero_state%valid_sort, size(aero_state%aero_weight), bin_grid, & 02403 all_procs_same) 02404 aero_state%valid_sort = .true. 02405 02406 end subroutine aero_state_sort 02407 02408 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02409 02410 !> Check that the sorted data is consistent. 02411 subroutine aero_state_check_sort(aero_state) 02412 02413 !> Aerosol state to check. 02414 type(aero_state_t), intent(in) :: aero_state 02415 02416 logical, parameter :: continue_on_error = .false. 02417 02418 integer :: i_part, i_bin 02419 02420 if (.not. aero_state%valid_sort) then 02421 write(0,*) 'SORTED CHECK ERROR: SORT NOT VALID' 02422 return 02423 end if 02424 02425 call aero_sorted_check(aero_state%aero_sorted, aero_state%apa, & 02426 size(aero_state%aero_weight), continue_on_error) 02427 02428 end subroutine aero_state_check_sort 02429 02430 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02431 02432 end module pmc_aero_state 02433