PartMC
2.1.5
|
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_bin_grid 00013 use pmc_aero_data 00014 use pmc_aero_particle 00015 use pmc_aero_dist 00016 use pmc_util 00017 use pmc_rand 00018 use pmc_aero_binned 00019 use pmc_mpi 00020 use pmc_spec_file 00021 use pmc_aero_info 00022 use pmc_aero_info_array 00023 use pmc_aero_weight 00024 #ifdef PMC_USE_MPI 00025 use mpi 00026 #endif 00027 00028 !> MPI tag for mixing particles between processes. 00029 integer, parameter :: AERO_STATE_TAG_MIX = 4987 00030 !> MPI tag for gathering between processes. 00031 integer, parameter :: AERO_STATE_TAG_GATHER = 4988 00032 !> MPI tag for scattering between processes. 00033 integer, parameter :: AERO_STATE_TAG_SCATTER = 4989 00034 00035 !> The current collection of aerosol particles. 00036 !! 00037 !! The particles in aero_state_t are stored sorted per-bin, to 00038 !! improve efficiency of access and sampling. If a particle has 00039 !! total radius \c r then calling <tt> i_bin = 00040 !! bin_grid_particle_in_bin(bin_grid, r)</tt> finds the bin number 00041 !! i_bin where that particle should go. That particle is then stored 00042 !! as \c aero_state%%bin(i_bin)%%particle(i_part), where \c i_part 00043 !! is the index within the bin. \c 00044 !! aero_state%%v(i_bin)%%p(i_part)%%vol(i_spec) is thus the volume 00045 !! of the \c i_spec-th species in the \c i_part-th particle in the 00046 !! \c i_bin-th bin. 00047 !! 00048 !! Typically most of the bins have only a few particles, while a 00049 !! small number of bins have many particles. To avoid having too 00050 !! much storage allocated for the bins with only a few particles, we 00051 !! do dynamic allocation and deallocation of the storage 00052 !! per-bin. With Fortran 90 we can't have arrays of arrays, so we 00053 !! have to use an array of pointers, and then allocate each pointer. 00054 !! 00055 !! To avoid doing allocation and deallocation every time we add or 00056 !! remove a particle to a bin, we always double or halve the bin 00057 !! storage as necessary. The actual number of particles stored in a 00058 !! bin will generally be less than the actual memory allocated for 00059 !! that bin, so we store the current number of particles in a bin in 00060 !! \c aero_state%%bin(i_bin)%%n_part. The allocated size of bin 00061 !! storage in \c aero_state%%bin(i_bin)%%particle is not stored 00062 !! explicitly, but can be obtained with the Fortran 90 SIZE() 00063 !! intrinsic function. 00064 !! 00065 !! Every time we remove particles we keep track of the particle ID 00066 !! and the action performed in the aero_info_array_t structure. This 00067 !! is typically cleared each time we output data to disk. 00068 type aero_state_t 00069 !> Bin arrays. 00070 type(aero_particle_array_t), pointer :: bin(:) 00071 !> Computational volume (m^3). 00072 real(kind=dp) :: comp_vol 00073 !> Total number of particles. 00074 integer :: n_part 00075 !> Information on removed particles. 00076 type(aero_info_array_t) :: aero_info_array 00077 end type aero_state_t 00078 00079 contains 00080 00081 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00082 00083 !> Allocates aerosol arrays. 00084 subroutine aero_state_allocate(aero_state) 00085 00086 !> Aerosol to initialize. 00087 type(aero_state_t), intent(out) :: aero_state 00088 00089 allocate(aero_state%bin(0)) 00090 call aero_info_array_allocate(aero_state%aero_info_array) 00091 00092 end subroutine aero_state_allocate 00093 00094 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00095 00096 !> Allocates aerosol arrays with the given sizes. 00097 subroutine aero_state_allocate_size(aero_state, n_bin, n_spec, n_source) 00098 00099 !> Aerosol to initialize. 00100 type(aero_state_t), intent(out) :: aero_state 00101 !> Number of bins. 00102 integer, intent(in) :: n_bin 00103 !> Number of species. 00104 integer, intent(in) :: n_spec 00105 !> Number of sources. 00106 integer, intent(in) :: n_source 00107 00108 integer i 00109 00110 allocate(aero_state%bin(n_bin)) 00111 do i = 1,n_bin 00112 call aero_particle_array_allocate_size(aero_state%bin(i), & 00113 0, n_spec, n_source) 00114 end do 00115 aero_state%comp_vol = 0d0 00116 aero_state%n_part = 0 00117 call aero_info_array_allocate(aero_state%aero_info_array) 00118 00119 end subroutine aero_state_allocate_size 00120 00121 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00122 00123 !> Deallocates a previously allocated aerosol. 00124 subroutine aero_state_deallocate(aero_state) 00125 00126 !> Aerosol to deallocate. 00127 type(aero_state_t), intent(inout) :: aero_state 00128 00129 integer :: n_bin, i 00130 00131 n_bin = size(aero_state%bin) 00132 do i = 1,n_bin 00133 call aero_particle_array_deallocate(aero_state%bin(i)) 00134 end do 00135 deallocate(aero_state%bin) 00136 call aero_info_array_deallocate(aero_state%aero_info_array) 00137 00138 end subroutine aero_state_deallocate 00139 00140 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00141 00142 !> Copies aerosol to a destination that has already had 00143 !> aero_state_allocate() called on it. 00144 subroutine aero_state_copy(aero_state_from, aero_state_to) 00145 00146 !> Reference aerosol. 00147 type(aero_state_t), intent(in) :: aero_state_from 00148 !> Already allocated. 00149 type(aero_state_t), intent(inout) :: aero_state_to 00150 00151 integer :: n_bin, i 00152 00153 n_bin = size(aero_state_from%bin) 00154 00155 call aero_state_deallocate(aero_state_to) 00156 call aero_state_allocate_size(aero_state_to, n_bin, 0, 0) 00157 00158 do i = 1,n_bin 00159 call aero_particle_array_copy(aero_state_from%bin(i), & 00160 aero_state_to%bin(i)) 00161 end do 00162 00163 aero_state_to%comp_vol = aero_state_from%comp_vol 00164 aero_state_to%n_part = aero_state_from%n_part 00165 00166 call aero_info_array_copy(aero_state_from%aero_info_array, & 00167 aero_state_to%aero_info_array) 00168 00169 end subroutine aero_state_copy 00170 00171 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00172 00173 !> Returns the total number of particles in an aerosol distribution. 00174 integer function aero_state_total_particles(aero_state) 00175 00176 !> Aerosol state. 00177 type(aero_state_t), intent(in) :: aero_state 00178 00179 aero_state_total_particles = aero_state%n_part 00180 00181 end function aero_state_total_particles 00182 00183 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00184 00185 !> Returns the total number of particles across all processes. 00186 integer function aero_state_total_particles_all_procs(aero_state) 00187 00188 !> Aerosol state. 00189 type(aero_state_t), intent(in) :: aero_state 00190 00191 call pmc_mpi_allreduce_sum_integer(& 00192 aero_state_total_particles(aero_state), & 00193 aero_state_total_particles_all_procs) 00194 00195 end function aero_state_total_particles_all_procs 00196 00197 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00198 00199 !> Resets an aero_state to have zero particles per bin. This must 00200 !> already have had aero_state_allocate() called on it. This 00201 !> function can be called more than once on the same state. 00202 subroutine aero_state_zero(aero_state) 00203 00204 !> State to zero. 00205 type(aero_state_t), intent(inout) :: aero_state 00206 00207 integer :: i, n_bin 00208 00209 n_bin = size(aero_state%bin) 00210 do i = 1,n_bin 00211 call aero_particle_array_zero(aero_state%bin(i)) 00212 end do 00213 aero_state%n_part = 0 00214 call aero_info_array_zero(aero_state%aero_info_array) 00215 00216 end subroutine aero_state_zero 00217 00218 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00219 00220 !> Add the given particle. 00221 subroutine aero_state_add_particle(aero_state, i_bin, aero_particle) 00222 00223 !> Aerosol state. 00224 type(aero_state_t), intent(inout) :: aero_state 00225 !> Bin number of particle to add. 00226 integer, intent(in) :: i_bin 00227 !> Particle to add. 00228 type(aero_particle_t), intent(in) :: aero_particle 00229 00230 call aero_particle_array_add_particle(aero_state%bin(i_bin), & 00231 aero_particle) 00232 aero_state%n_part = aero_state%n_part + 1 00233 00234 end subroutine aero_state_add_particle 00235 00236 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00237 00238 !> Remove the given particle without recording it. 00239 subroutine aero_state_remove_particle_no_info(aero_state, i_bin, index) 00240 00241 !> Aerosol state. 00242 type(aero_state_t), intent(inout) :: aero_state 00243 !> Bin number of particle to remove. 00244 integer, intent(in) :: i_bin 00245 !> Index in bin of particle to remove. 00246 integer, intent(in) :: index 00247 00248 call aero_particle_array_remove_particle(aero_state%bin(i_bin), index) 00249 aero_state%n_part = aero_state%n_part - 1 00250 00251 end subroutine aero_state_remove_particle_no_info 00252 00253 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00254 00255 !> Remove the given particle and record the removal. 00256 subroutine aero_state_remove_particle_with_info(aero_state, i_bin, & 00257 index, aero_info) 00258 00259 !> Aerosol state. 00260 type(aero_state_t), intent(inout) :: aero_state 00261 !> Bin number of particle to remove. 00262 integer, intent(in) :: i_bin 00263 !> Index in bin of particle to remove. 00264 integer, intent(in) :: index 00265 !> Removal info. 00266 type(aero_info_t), intent(in) :: aero_info 00267 00268 call aero_state_remove_particle_no_info(aero_state, i_bin, index) 00269 call aero_info_array_add_aero_info(aero_state%aero_info_array, & 00270 aero_info) 00271 00272 end subroutine aero_state_remove_particle_with_info 00273 00274 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00275 00276 !> Remove the given particle and possibly record the removal. 00277 subroutine aero_state_remove_particle(aero_state, i_bin, index, & 00278 record_removal, aero_info) 00279 00280 !> Aerosol state. 00281 type(aero_state_t), intent(inout) :: aero_state 00282 !> Bin number of particle to remove. 00283 integer, intent(in) :: i_bin 00284 !> Index in bin of particle to remove. 00285 integer, intent(in) :: index 00286 !> Whether to record the removal in the aero_info_array. 00287 logical, intent(in) :: record_removal 00288 !> Removal info. 00289 type(aero_info_t), intent(in) :: aero_info 00290 00291 if (record_removal) then 00292 call aero_state_remove_particle_with_info(aero_state, i_bin, & 00293 index, aero_info) 00294 else 00295 call aero_state_remove_particle_no_info(aero_state, i_bin, index) 00296 end if 00297 00298 end subroutine aero_state_remove_particle 00299 00300 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00301 00302 subroutine aero_state_remove_rand_particle_from_bin(aero_state, & 00303 i_bin, aero_particle) 00304 00305 !> Aerosol state. 00306 type(aero_state_t), intent(inout) :: aero_state 00307 !> Bin number to remove particle from. 00308 integer, intent(in) :: i_bin 00309 !> Removed particle. 00310 type(aero_particle_t), intent(inout) :: aero_particle 00311 00312 integer :: i_part 00313 00314 call assert(392182617, aero_state%bin(i_bin)%n_part > 0) 00315 i_part = pmc_rand_int(aero_state%bin(i_bin)%n_part) 00316 call aero_particle_copy(aero_state%bin(i_bin)%particle(i_part), & 00317 aero_particle) 00318 call aero_state_remove_particle_no_info(aero_state, i_bin, i_part) 00319 00320 end subroutine aero_state_remove_rand_particle_from_bin 00321 00322 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00323 00324 !> <tt>aero_state += aero_state_delta</tt>, with adding the computational 00325 !> volumes, so the new concentration is the (volume-weighted) 00326 !> average of the two concentration. 00327 subroutine aero_state_add(aero_state, aero_state_delta) 00328 00329 !> Aerosol state. 00330 type(aero_state_t), intent(inout) :: aero_state 00331 !> Increment. 00332 type(aero_state_t), intent(in) :: aero_state_delta 00333 00334 call aero_state_add_particles(aero_state, aero_state_delta) 00335 aero_state%comp_vol = aero_state%comp_vol + aero_state_delta%comp_vol 00336 00337 end subroutine aero_state_add 00338 00339 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00340 00341 !> <tt>aero_state += aero_state_delta</tt>, with the computational 00342 !> volume of \c aero_state left unchanged, so the new concentration is the 00343 !> sum of the two concentrations, computed with \c aero_state%comp_vol. 00344 subroutine aero_state_add_particles(aero_state, aero_state_delta) 00345 00346 !> Aerosol state. 00347 type(aero_state_t), intent(inout) :: aero_state 00348 !> Increment. 00349 type(aero_state_t), intent(in) :: aero_state_delta 00350 00351 integer :: i_bin, i_part 00352 00353 call assert(265083067, & 00354 size(aero_state%bin) == size(aero_state_delta%bin)) 00355 do i_bin = 1,size(aero_state_delta%bin) 00356 do i_part = 1,aero_state_delta%bin(i_bin)%n_part 00357 call aero_state_add_particle(aero_state, i_bin, & 00358 aero_state_delta%bin(i_bin)%particle(i_part)) 00359 end do 00360 end do 00361 call aero_info_array_add(aero_state%aero_info_array, & 00362 aero_state_delta%aero_info_array) 00363 00364 end subroutine aero_state_add_particles 00365 00366 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00367 00368 !> Generates a Poisson sample of an aero_dist, adding to 00369 !> aero_state. The sampled amount is sample_prop * 00370 !> aero_state%comp_vol. 00371 subroutine aero_state_add_aero_dist_sample(aero_state, bin_grid, & 00372 aero_data, aero_weight, aero_dist, sample_prop, create_time) 00373 00374 !> Aero state to add to. 00375 type(aero_state_t), intent(inout) :: aero_state 00376 !> Bin grid. 00377 type(bin_grid_t), intent(in) :: bin_grid 00378 !> Aero data values. 00379 type(aero_data_t), intent(in) :: aero_data 00380 !> Aerosol weight. 00381 type(aero_weight_t), intent(in) :: aero_weight 00382 !> Distribution to sample. 00383 type(aero_dist_t), intent(in) :: aero_dist 00384 !> Volume fraction to sample (1). 00385 real(kind=dp), intent(in) :: sample_prop 00386 !> Creation time for new particles (s). 00387 real(kind=dp), intent(in) :: create_time 00388 00389 real(kind=dp) :: n_samp_avg, sample_vol, radius, vol 00390 integer :: n_samp, i_mode, i_samp, i_bin 00391 integer :: num_per_bin(bin_grid%n_bin) 00392 type(aero_mode_t), pointer :: aero_mode 00393 type(aero_particle_t) :: aero_particle 00394 00395 call aero_particle_allocate_size(aero_particle, aero_data%n_spec, & 00396 aero_data%n_source) 00397 sample_vol = sample_prop * aero_state%comp_vol 00398 do i_mode = 1,aero_dist%n_mode 00399 aero_mode => aero_dist%mode(i_mode) 00400 n_samp_avg = sample_vol & 00401 * aero_mode_weighted_num_conc(aero_mode, aero_weight) 00402 n_samp = rand_poisson(n_samp_avg) 00403 do i_samp = 1,n_samp 00404 call aero_particle_zero(aero_particle) 00405 call aero_mode_sample_radius(aero_mode, aero_weight, radius) 00406 vol = rad2vol(radius) 00407 call aero_particle_set_vols(aero_particle, aero_mode%vol_frac * vol) 00408 call aero_particle_new_id(aero_particle) 00409 call aero_particle_set_create_time(aero_particle, create_time) 00410 call aero_particle_set_source(aero_particle, aero_mode%source) 00411 i_bin = aero_particle_in_bin(aero_particle, bin_grid) 00412 call aero_state_add_particle(aero_state, i_bin, aero_particle) 00413 end do 00414 end do 00415 call aero_particle_deallocate(aero_particle) 00416 00417 end subroutine aero_state_add_aero_dist_sample 00418 00419 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00420 00421 !> Choose a random particle from the aero_state. 00422 subroutine aero_state_rand_particle(aero_state, i_bin, i_part) 00423 00424 !> Original state. 00425 type(aero_state_t), intent(in) :: aero_state 00426 !> Bin number of particle. 00427 integer, intent(out) :: i_bin 00428 !> Particle number within bin. 00429 integer, intent(out) :: i_part 00430 00431 integer :: n_bin, disc_pdf(size(aero_state%bin)) 00432 00433 call assert(950725003, aero_state%n_part > 0) 00434 n_bin = size(aero_state%bin) 00435 disc_pdf = (/(aero_state%bin(i_bin)%n_part, i_bin = 1,n_bin)/) 00436 i_bin = sample_disc_pdf(disc_pdf) 00437 i_part = pmc_rand_int(aero_state%bin(i_bin)%n_part) 00438 00439 end subroutine aero_state_rand_particle 00440 00441 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00442 00443 !> Generates a random sample by removing particles from 00444 !> aero_state_from and adding them to aero_state_to, which must 00445 !> be already allocated (and should have its comp_vol set). 00446 subroutine aero_state_sample(aero_state_from, aero_state_to, & 00447 sample_prob, removal_action) 00448 00449 !> Original state. 00450 type(aero_state_t), intent(inout) :: aero_state_from 00451 !> Destination state. 00452 type(aero_state_t), intent(inout) :: aero_state_to 00453 !> Probability of sampling each particle. 00454 real(kind=dp), intent(in) :: sample_prob 00455 !> Action for removal (see pmc_aero_info module for action 00456 !> parameters). Set to AERO_INFO_NONE to not log removal. 00457 integer, intent(in) :: removal_action 00458 00459 integer :: n_transfer, i_transfer, n_bin, i_bin, i_part 00460 logical :: do_add, do_remove 00461 real(kind=dp) :: vol_ratio 00462 type(aero_info_t) :: aero_info 00463 00464 call assert(721006962, (sample_prob >= 0d0) .and. (sample_prob <= 1d0)) 00465 n_transfer = rand_binomial(aero_state_total_particles(aero_state_from), & 00466 sample_prob) 00467 n_bin = size(aero_state_from%bin) 00468 vol_ratio = aero_state_to%comp_vol / aero_state_from%comp_vol 00469 i_transfer = 0 00470 do while (i_transfer < n_transfer) 00471 if (aero_state_total_particles(aero_state_from) <= 0) exit 00472 call aero_state_rand_particle(aero_state_from, i_bin, i_part) 00473 if (aero_state_to%comp_vol == aero_state_from%comp_vol) then 00474 ! to_comp_vol == from_comp_vol so just move the particle 00475 do_add = .true. 00476 do_remove = .true. 00477 elseif (aero_state_to%comp_vol > aero_state_from%comp_vol) then 00478 ! to_comp_vol is bigger than from_comp_vol, so only maybe 00479 ! remove the particle but always add it (vol_ratio > 1) 00480 do_add = .true. 00481 do_remove = .false. 00482 if (pmc_random() < 1d0 / vol_ratio) then 00483 do_remove = .true. 00484 end if 00485 else ! aero_state_to%comp_vol < aero_state_from%comp_vol 00486 ! to_comp_vol is smaller than from_comp_vol, so always 00487 ! remove the particle but only maybe add it (vol_ratio < 1) 00488 do_add = .false. 00489 if (pmc_random() < vol_ratio) then 00490 do_add = .true. 00491 end if 00492 do_remove = .true. 00493 end if 00494 if (do_add) then 00495 call aero_state_add_particle(aero_state_to, i_bin, & 00496 aero_state_from%bin(i_bin)%particle(i_part)) 00497 end if 00498 if (do_remove) then 00499 if (removal_action /= AERO_INFO_NONE) then 00500 call aero_info_allocate(aero_info) 00501 aero_info%id = & 00502 aero_state_from%bin(i_bin)%particle(i_part)%id 00503 aero_info%action = removal_action 00504 call aero_state_remove_particle(aero_state_from, i_bin, & 00505 i_part, .true., aero_info) 00506 call aero_info_deallocate(aero_info) 00507 else 00508 call aero_state_remove_particle(aero_state_from, i_bin, & 00509 i_part, .false., aero_info) 00510 end if 00511 i_transfer = i_transfer + 1 00512 end if 00513 end do 00514 00515 end subroutine aero_state_sample 00516 00517 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00518 00519 !> Create the bin number and mass arrays from aero_state%v. 00520 subroutine aero_state_to_binned(bin_grid, aero_data, aero_weight, & 00521 aero_state, aero_binned) 00522 00523 !> Bin grid. 00524 type(bin_grid_t), intent(in) :: bin_grid 00525 !> Aerosol data. 00526 type(aero_data_t), intent(in) :: aero_data 00527 !> Aerosol weight. 00528 type(aero_weight_t), intent(in) :: aero_weight 00529 !> Aerosol state. 00530 type(aero_state_t), intent(in) :: aero_state 00531 !> Binned distributions. 00532 type(aero_binned_t), intent(inout) :: aero_binned 00533 00534 integer :: b, j, s 00535 type(aero_particle_t), pointer :: aero_particle 00536 00537 aero_binned%num_conc = 0d0 00538 aero_binned%vol_conc = 0d0 00539 do b = 1,bin_grid%n_bin 00540 do j = 1,aero_state%bin(b)%n_part 00541 aero_particle => aero_state%bin(b)%particle(j) 00542 aero_binned%vol_conc(b,:) = aero_binned%vol_conc(b,:) & 00543 + aero_particle%vol / aero_state%comp_vol & 00544 * aero_weight_value(aero_weight, & 00545 aero_particle_radius(aero_particle)) / bin_grid%log_width 00546 aero_binned%num_conc(b) = aero_binned%num_conc(b) & 00547 + 1d0 / aero_state%comp_vol & 00548 * aero_weight_value(aero_weight, & 00549 aero_particle_radius(aero_particle)) / bin_grid%log_width 00550 end do 00551 end do 00552 00553 end subroutine aero_state_to_binned 00554 00555 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00556 00557 !> Does the same thing as aero_state_to_bin() but based on dry radius. 00558 subroutine aero_state_to_binned_dry(bin_grid, aero_data, aero_weight, & 00559 aero_state, aero_binned) 00560 00561 !> Bin grid. 00562 type(bin_grid_t), intent(in) :: bin_grid 00563 !> Aerosol data. 00564 type(aero_data_t), intent(in) :: aero_data 00565 !> Aerosol weight. 00566 type(aero_weight_t), intent(in) :: aero_weight 00567 !> Aerosol state. 00568 type(aero_state_t), intent(in) :: aero_state 00569 !> Binned distributions. 00570 type(aero_binned_t), intent(inout) :: aero_binned 00571 00572 integer :: b, j, s, b_dry 00573 type(aero_particle_t), pointer :: aero_particle 00574 00575 aero_binned%num_conc = 0d0 00576 aero_binned%vol_conc = 0d0 00577 do b = 1,bin_grid%n_bin 00578 do j = 1,aero_state%bin(b)%n_part 00579 aero_particle => aero_state%bin(b)%particle(j) 00580 b_dry = bin_grid_particle_in_bin(bin_grid, & 00581 aero_particle_solute_radius(aero_particle, aero_data)) 00582 aero_binned%vol_conc(b_dry,:) = aero_binned%vol_conc(b_dry,:) & 00583 + aero_particle%vol / aero_state%comp_vol & 00584 * aero_weight_value(aero_weight, & 00585 aero_particle_radius(aero_particle)) / bin_grid%log_width 00586 aero_binned%num_conc(b_dry) = aero_binned%num_conc(b_dry) & 00587 + 1d0 / aero_state%comp_vol & 00588 * aero_weight_value(aero_weight, & 00589 aero_particle_radius(aero_particle)) / bin_grid%log_width 00590 end do 00591 end do 00592 00593 end subroutine aero_state_to_binned_dry 00594 00595 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00596 00597 !> Doubles number of particles. 00598 subroutine aero_state_double(aero_state) 00599 00600 !> Aerosol state. 00601 type(aero_state_t), intent(inout) :: aero_state 00602 00603 integer :: i, n_bin 00604 00605 n_bin = size(aero_state%bin) 00606 do i = 1,n_bin 00607 call aero_particle_array_double(aero_state%bin(i)) 00608 end do 00609 aero_state%comp_vol = 2d0 * aero_state%comp_vol 00610 aero_state%n_part = 2 * aero_state%n_part 00611 00612 end subroutine aero_state_double 00613 00614 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00615 00616 !> Remove approximately half of the particles in each bin. 00617 subroutine aero_state_halve(aero_state, bin_grid) 00618 00619 !> Aerosol state. 00620 type(aero_state_t), intent(inout) :: aero_state 00621 !> Bin grid. 00622 type(bin_grid_t), intent(in) :: bin_grid 00623 00624 integer :: i_bin, i_part, n_part_orig, i_remove, n_remove 00625 type(aero_info_t) :: aero_info 00626 00627 n_part_orig = aero_state%n_part 00628 do i_bin = 1,bin_grid%n_bin 00629 n_remove = prob_round(real(aero_state%bin(i_bin)%n_part, kind=dp) 00630 / 2d0) 00631 do i_remove = 1,n_remove 00632 i_part = pmc_rand_int(aero_state%bin(i_bin)%n_part) 00633 call aero_info_allocate(aero_info) 00634 aero_info%id = & 00635 aero_state%bin(i_bin)%particle(i_part)%id 00636 aero_info%action = AERO_INFO_HALVED 00637 call aero_state_remove_particle(aero_state, i_bin, & 00638 i_part, .true., aero_info) 00639 call aero_info_deallocate(aero_info) 00640 end do 00641 end do 00642 aero_state%comp_vol = aero_state%comp_vol & 00643 * real(aero_state%n_part, kind=dp) / real(n_part_orig, kind=dp) 00644 00645 end subroutine aero_state_halve 00646 00647 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00648 00649 !> Takes an aero_state_t where the particles might no longer be in 00650 !> the correct bins and resorts it so that every particle is in the 00651 !> correct bin. 00652 subroutine aero_state_resort(bin_grid, aero_state) 00653 00654 !> Bin_grid. 00655 type(bin_grid_t), intent(in) :: bin_grid 00656 !> Aerosol state. 00657 type(aero_state_t), intent(inout) :: aero_state 00658 00659 integer :: i_bin, j, i_new_bin, k 00660 type(aero_info_t) :: aero_info ! dummy variable, never used 00661 00662 ! The approach here is inefficient because we might reprocess 00663 ! particles. For example, if we are doing bin 1 and we shift a 00664 ! particle up to bin 2, when we do bin 2 we will reprocess it. It 00665 ! seems to be more trouble than it's worth to worry about this, 00666 ! however. 00667 00668 do i_bin = 1,bin_grid%n_bin 00669 j = 1 00670 do while (j .le. aero_state%bin(i_bin)%n_part) 00671 ! find the new bin 00672 i_new_bin = aero_particle_in_bin( & 00673 aero_state%bin(i_bin)%particle(j), bin_grid) 00674 00675 ! if the bin number has changed, move the particle 00676 if (i_bin .ne. i_new_bin) then 00677 call aero_state_add_particle(aero_state, i_new_bin, & 00678 aero_state%bin(i_bin)%particle(j)) 00679 call aero_state_remove_particle(aero_state, i_bin, j, & 00680 .false., aero_info) 00681 00682 ! in this case, don't advance j, so that we will still 00683 ! process the particle we just moved into the hole 00684 else 00685 ! if we didn't move the particle, advance j to process 00686 ! the next particle 00687 j = j + 1 00688 end if 00689 end do 00690 end do 00691 00692 ! now shrink the bin storage if necessary 00693 do i_bin = 1,bin_grid%n_bin 00694 call aero_particle_array_shrink(aero_state%bin(i_bin)) 00695 end do 00696 00697 end subroutine aero_state_resort 00698 00699 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00700 00701 !> Mix the aero_states between all processes. Currently uses a 00702 !> simple all-to-all diffusion. 00703 subroutine aero_state_mix(aero_state, del_t, mix_timescale, & 00704 aero_data, bin_grid) 00705 00706 !> Aerosol state. 00707 type(aero_state_t), intent(inout) :: aero_state 00708 !> Timestep (s). 00709 real(kind=dp), intent(in) :: del_t 00710 !> Mixing timescale (s). 00711 real(kind=dp), intent(in) :: mix_timescale 00712 !> Aero data values. 00713 type(aero_data_t), intent(in) :: aero_data 00714 !> Bin grid. 00715 type(bin_grid_t), intent(in) :: bin_grid 00716 00717 #ifdef PMC_USE_MPI 00718 integer :: rank, n_proc, dest_proc, ierr 00719 integer :: buffer_size, buffer_size_check 00720 character, allocatable :: buffer(:) 00721 type(aero_state_t), allocatable :: aero_state_sends(:) 00722 real(kind=dp) :: prob_transfer, prob_not_transferred 00723 real(kind=dp) :: prob_transfer_given_not_transferred 00724 real(kind=dp), allocatable :: comp_vols(:) 00725 00726 ! process information 00727 rank = pmc_mpi_rank() 00728 n_proc = pmc_mpi_size() 00729 if (n_proc == 1) then 00730 ! buffer allocation below fails if n_proc == 1 00731 ! so bail out early (nothing to mix anyway) 00732 return 00733 end if 00734 allocate(comp_vols(n_proc)) 00735 call mpi_allgather(aero_state%comp_vol, 1, MPI_REAL8, & 00736 comp_vols, 1, MPI_REAL8, MPI_COMM_WORLD, ierr) 00737 00738 ! extract particles to send 00739 allocate(aero_state_sends(n_proc)) 00740 prob_not_transferred = 1d0 00741 do dest_proc = 0,(n_proc - 1) 00742 if (dest_proc /= rank) then 00743 call aero_state_allocate_size(aero_state_sends(dest_proc + 1), & 00744 bin_grid%n_bin, aero_data%n_spec, aero_data%n_source) 00745 aero_state_sends(dest_proc + 1)%comp_vol = aero_state%comp_vol 00746 prob_transfer = (1d0 - exp(- del_t / mix_timescale)) & 00747 / real(n_proc, kind=dp) 00748 * min(1d0, comp_vols(dest_proc + 1) / comp_vols(rank + 1)) 00749 ! because we are doing sequential sampling from the aero_state 00750 ! we need to scale up the later transfer probabilities, because 00751 ! the later particles are being transferred conditioned on the 00752 ! fact that they did not transfer already 00753 prob_transfer_given_not_transferred = prob_transfer & 00754 / prob_not_transferred 00755 call aero_state_sample(aero_state, & 00756 aero_state_sends(dest_proc + 1), & 00757 prob_transfer_given_not_transferred, AERO_INFO_NONE) 00758 prob_not_transferred = prob_not_transferred - prob_transfer 00759 end if 00760 end do 00761 00762 ! send the particles 00763 buffer_size = 0 00764 do dest_proc = 0,(n_proc - 1) 00765 if (dest_proc /= rank) then 00766 buffer_size = buffer_size & 00767 + pmc_mpi_pack_size_aero_state(aero_state_sends(dest_proc + 1)) 00768 buffer_size = buffer_size + MPI_BSEND_OVERHEAD + 1024 00769 end if 00770 end do 00771 allocate(buffer(buffer_size)) 00772 call mpi_buffer_attach(buffer, buffer_size, ierr) 00773 call pmc_mpi_check_ierr(ierr) 00774 do dest_proc = 0,(n_proc - 1) 00775 if (dest_proc /= rank) then 00776 call send_aero_state_mix(dest_proc, aero_state_sends(dest_proc + 1)) 00777 call aero_state_deallocate(aero_state_sends(dest_proc + 1)) 00778 end if 00779 end do 00780 deallocate(aero_state_sends) 00781 00782 ! receive the particles 00783 do dest_proc = 0,(n_proc - 1) 00784 if (dest_proc /= rank) then 00785 call recv_aero_state_mix(aero_state) 00786 end if 00787 end do 00788 00789 ! clean up 00790 call pmc_mpi_barrier() ! syncronize to make sure all buffers are done 00791 call mpi_buffer_detach(buffer, buffer_size_check, ierr) 00792 call pmc_mpi_check_ierr(ierr) 00793 call assert(995066222, buffer_size == buffer_size_check) 00794 deallocate(comp_vols) 00795 #endif 00796 00797 end subroutine aero_state_mix 00798 00799 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00800 00801 !> Send the given \c aero_state to the destination process for 00802 !> mixing. 00803 subroutine send_aero_state_mix(dest_proc, aero_state) 00804 00805 !> Destination process number. 00806 integer, intent(in) :: dest_proc 00807 !> Aero_state to send. 00808 type(aero_state_t), intent(in) :: aero_state 00809 00810 #ifdef PMC_USE_MPI 00811 character, allocatable :: buffer(:) 00812 integer :: buffer_size, max_buffer_size, position, ierr 00813 00814 max_buffer_size = 0 00815 max_buffer_size = max_buffer_size & 00816 + pmc_mpi_pack_size_aero_state(aero_state) 00817 allocate(buffer(max_buffer_size)) 00818 position = 0 00819 call pmc_mpi_pack_aero_state(buffer, position, aero_state) 00820 call assert(311031745, position <= max_buffer_size) 00821 buffer_size = position 00822 call mpi_bsend(buffer, buffer_size, MPI_CHARACTER, dest_proc, & 00823 AERO_STATE_TAG_MIX, MPI_COMM_WORLD, ierr) 00824 call pmc_mpi_check_ierr(ierr) 00825 deallocate(buffer) 00826 #endif 00827 00828 end subroutine send_aero_state_mix 00829 00830 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00831 00832 !> Receive exactly one \c aero_state for mixing and add it on to the 00833 !> given \c aero_state. 00834 subroutine recv_aero_state_mix(aero_state) 00835 00836 !> Base \c aero_state to add new particles to. 00837 type(aero_state_t), intent(inout) :: aero_state 00838 00839 #ifdef PMC_USE_MPI 00840 integer :: buffer_size, check_buffer_size, position, sent_proc, ierr 00841 integer :: status(MPI_STATUS_SIZE) 00842 character, allocatable :: buffer(:) 00843 type(aero_state_t) :: aero_state_mix 00844 00845 ! get the message size 00846 call mpi_probe(MPI_ANY_SOURCE, AERO_STATE_TAG_MIX, MPI_COMM_WORLD, & 00847 status, ierr) 00848 call pmc_mpi_check_ierr(ierr) 00849 sent_proc = status(MPI_SOURCE) 00850 call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr) 00851 call pmc_mpi_check_ierr(ierr) 00852 allocate(buffer(buffer_size)) 00853 00854 ! get the message 00855 call mpi_recv(buffer, buffer_size, MPI_CHARACTER, & 00856 sent_proc, AERO_STATE_TAG_MIX, MPI_COMM_WORLD, status, ierr) 00857 call pmc_mpi_check_ierr(ierr) 00858 call mpi_get_count(status, MPI_CHARACTER, check_buffer_size, ierr) 00859 call pmc_mpi_check_ierr(ierr) 00860 call assert(967116476, buffer_size == check_buffer_size) 00861 call assert(377784810, sent_proc == status(MPI_SOURCE)) 00862 00863 ! unpack it 00864 position = 0 00865 call aero_state_allocate(aero_state_mix) 00866 call pmc_mpi_unpack_aero_state(buffer, position, aero_state_mix) 00867 call assert(908434410, position == buffer_size) 00868 deallocate(buffer) 00869 00870 ! add the particles 00871 call aero_state_add_particles(aero_state, aero_state_mix) 00872 call aero_state_deallocate(aero_state_mix) 00873 #endif 00874 00875 end subroutine recv_aero_state_mix 00876 00877 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00878 00879 !> Check that all particles are in the correct bins and that the 00880 !> bin numbers and masses are correct. This is for debugging only. 00881 subroutine aero_state_check(bin_grid, aero_data, aero_state) 00882 00883 !> Bin_grid. 00884 type(bin_grid_t), intent(in) :: bin_grid 00885 !> Aerosol data. 00886 type(aero_data_t), intent(in) :: aero_data 00887 !> Aerosol state. 00888 type(aero_state_t), intent(inout) :: aero_state 00889 00890 real(kind=dp) :: check_bin_v, check_vol_conc, vol_tol 00891 real(kind=dp) :: num_tol, state_num_conc 00892 integer :: i, k, k_check, s, n_part_check, id, max_id 00893 logical :: error 00894 logical, allocatable :: id_present(:) 00895 00896 error = .false. 00897 00898 ! check that the total number of particles is correct 00899 i = 0 ! HACK to shut up gfortran warning 00900 n_part_check = sum((/(aero_state%bin(i)%n_part, i = 1,bin_grid%n_bin)/)) 00901 if (aero_state%n_part /= n_part_check) then 00902 write(0,'(a20,a20)') 'aero_state%n_part', 'n_part_check' 00903 write(0,'(i20,i20)') aero_state%n_part, n_part_check 00904 error = .true. 00905 end if 00906 00907 ! check that all particles are in the correct bins 00908 do k = 1,bin_grid%n_bin 00909 do i = 1,aero_state%bin(k)%n_part 00910 k_check = aero_particle_in_bin(aero_state%bin(k)%particle(i), & 00911 bin_grid) 00912 if (k .ne. k_check) then 00913 write(0,'(a10,a10,a10)') 'i', 'k', 'k_check' 00914 write(0,'(i10,i10,i10)') i, k, k_check 00915 error = .true. 00916 end if 00917 end do 00918 end do 00919 00920 ! check we don't have duplicate IDs 00921 max_id = 0 00922 do k = 1,bin_grid%n_bin 00923 do i = 1,aero_state%bin(k)%n_part 00924 id = aero_state%bin(k)%particle(i)%id 00925 if (id > max_id) max_id = id 00926 end do 00927 end do 00928 allocate(id_present(max_id)) 00929 do i = 1,max_id 00930 id_present(i) = .false. 00931 end do 00932 do k = 1,bin_grid%n_bin 00933 do i = 1,aero_state%bin(k)%n_part 00934 id = aero_state%bin(k)%particle(i)%id 00935 if (id_present(id)) then 00936 write(0,'(a15,a10,a10)') 'duplicate id', 'bin', 'index' 00937 write(0,'(i15,i10,i10)') id, k, i 00938 error = .true. 00939 end if 00940 id_present(id) = .true. 00941 end do 00942 end do 00943 deallocate(id_present) 00944 00945 if (error) then 00946 call die_msg(371923719, 'aero_state_check() failed') 00947 end if 00948 00949 end subroutine aero_state_check 00950 00951 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00952 00953 !> Set each aerosol particle to have its original total volume, but 00954 !> species volume ratios given by the total species volume ratio 00955 !> within each bin. This preserves the (weighted) total species 00956 !> volume per bin as well as per-particle total volumes. 00957 subroutine aero_state_bin_average_comp(aero_state, bin_grid, aero_data, & 00958 aero_weight, dry_volume) 00959 00960 !> Aerosol state to average. 00961 type(aero_state_t), intent(inout) :: aero_state 00962 !> Bin grid to average within. 00963 type(bin_grid_t), intent(in) :: bin_grid 00964 !> Aerosol data. 00965 type(aero_data_t), intent(in) :: aero_data 00966 !> Aerosol weight. 00967 type(aero_weight_t), intent(in) :: aero_weight 00968 !> Whether to use dry volume (rather than wet). 00969 logical, intent(in) :: dry_volume 00970 00971 real(kind=dp) :: weighted_species_volumes(aero_data%n_spec) 00972 real(kind=dp) :: weighted_total_volume, particle_volume, weight 00973 integer :: i_bin, i_part, i_spec 00974 type(aero_particle_t), pointer :: aero_particle 00975 00976 do i_bin = 1,bin_grid%n_bin 00977 weighted_species_volumes = 0d0 00978 weighted_total_volume = 0d0 00979 do i_part = 1,aero_state%bin(i_bin)%n_part 00980 aero_particle => aero_state%bin(i_bin)%particle(i_part) 00981 weight = aero_weight_value(aero_weight, & 00982 aero_particle_radius(aero_particle)) 00983 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 00984 aero_data, dry_volume) 00985 weighted_species_volumes = weighted_species_volumes & 00986 + weight * aero_particle%vol 00987 weighted_total_volume = weighted_total_volume & 00988 + weight * particle_volume 00989 end do 00990 do i_part = 1,aero_state%bin(i_bin)%n_part 00991 aero_particle => aero_state%bin(i_bin)%particle(i_part) 00992 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 00993 aero_data, dry_volume) 00994 aero_particle%vol = particle_volume * weighted_species_volumes & 00995 / weighted_total_volume 00996 if (dry_volume .and. (aero_data%i_water > 0)) then 00997 ! set water to zero if we are doing dry volume averaging 00998 aero_particle%vol(aero_data%i_water) = 0d0 00999 end if 01000 end do 01001 end do 01002 01003 end subroutine aero_state_bin_average_comp 01004 01005 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01006 01007 !> Set each aerosol particle to have its original species ratios, 01008 !> but total volume given by the average volume of all particles 01009 !> within each bin. 01010 !! 01011 !! This does not preserve the total species volume 01012 !! per bin. If the \c bin_center parameter is \c .true. then the 01013 !! particles in each bin are set to have the bin center volume, 01014 !! rather than the average volume of the particles in that bin. 01015 !! 01016 !! If the weighting function is not constant (AERO_WEIGHT_TYPE_NONE) 01017 !! then the averaging can be performed in either a number-preserving 01018 !! way or in a volume-preserving way. The volume-preserving way does 01019 !! not preserve species volume ratios in gernal, but will do so if 01020 !! the particle population has already been composition-averaged. 01021 subroutine aero_state_bin_average_size(aero_state, bin_grid, aero_data, & 01022 aero_weight, dry_volume, bin_center, preserve_number) 01023 01024 !> Aerosol state to average. 01025 type(aero_state_t), intent(inout) :: aero_state 01026 !> Bin grid to average within. 01027 type(bin_grid_t), intent(in) :: bin_grid 01028 !> Aerosol data. 01029 type(aero_data_t), intent(in) :: aero_data 01030 !> Aerosol weight. 01031 type(aero_weight_t), intent(in) :: aero_weight 01032 !> Whether to use dry volume (rather than wet). 01033 logical, intent(in) :: dry_volume 01034 !> Whether to assign the bin center volume (rather than the average 01035 !> volume). 01036 logical, intent(in) :: bin_center 01037 !> Whether to use the number-preserving scheme (otherwise will use 01038 !> the volume-preserving scheme). This parameter has no effect if 01039 !> \c bin_center is \c .true. 01040 logical, intent(in) :: preserve_number 01041 01042 real(kind=dp) :: weighted_total_volume, particle_volume 01043 real(kind=dp) :: new_particle_volume, weight, total_weight 01044 real(kind=dp) :: lower_volume, upper_volume, center_volume 01045 real(kind=dp) :: lower_function, upper_function, center_function 01046 integer :: i_bin, i_part, i_bisect 01047 type(aero_particle_t), pointer :: aero_particle 01048 01049 do i_bin = 1,bin_grid%n_bin 01050 if (aero_state%bin(i_bin)%n_part == 0) then 01051 cycle 01052 end if 01053 01054 total_weight = 0d0 01055 weighted_total_volume = 0d0 01056 do i_part = 1,aero_state%bin(i_bin)%n_part 01057 aero_particle => aero_state%bin(i_bin)%particle(i_part) 01058 weight = aero_weight_value(aero_weight, & 01059 aero_particle_radius(aero_particle)) 01060 total_weight = total_weight + weight 01061 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01062 aero_data, dry_volume) 01063 weighted_total_volume = weighted_total_volume & 01064 + weight * particle_volume 01065 end do 01066 01067 ! determine the new_particle_volume for all particles in this bin 01068 if (bin_center) then 01069 new_particle_volume = rad2vol(bin_grid%center_radius(i_bin)) 01070 elseif (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then 01071 new_particle_volume = weighted_total_volume & 01072 / real(aero_state%bin(i_bin)%n_part, kind=dp) 01073 elseif (preserve_number) then 01074 ! number-preserving scheme: Solve the implicit equation: 01075 ! n_part * W(new_vol) = total_weight 01076 ! 01077 ! We assume that the weighting function is strictly monotone 01078 ! so this equation has a unique solution and the solution 01079 ! lies between the min and max particle volumes. We use 01080 ! bisection as this doesn't really need to be fast, just 01081 ! robust. 01082 01083 ! initialize to min and max particle volumes 01084 do i_part = 1,aero_state%bin(i_bin)%n_part 01085 aero_particle => aero_state%bin(i_bin)%particle(i_part) 01086 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01087 aero_data, dry_volume) 01088 if (i_part == 1) then 01089 lower_volume = particle_volume 01090 upper_volume = particle_volume 01091 end if 01092 lower_volume = min(lower_volume, particle_volume) 01093 upper_volume = max(upper_volume, particle_volume) 01094 end do 01095 01096 lower_function = real(aero_state%bin(i_bin)%n_part, kind=dp) 01097 * aero_weight_value(aero_weight, vol2rad(lower_volume)) 01098 - total_weight 01099 upper_function = real(aero_state%bin(i_bin)%n_part, kind=dp) 01100 * aero_weight_value(aero_weight, vol2rad(upper_volume)) 01101 - total_weight 01102 01103 ! do 50 rounds of bisection (2^50 = 10^15) 01104 do i_bisect = 1,50 01105 center_volume = (lower_volume + upper_volume) / 2d0 01106 center_function = real(aero_state%bin(i_bin)%n_part, kind=dp) 01107 * aero_weight_value(aero_weight, vol2rad(center_volume)) 01108 - total_weight 01109 if ((lower_function > 0d0 .and. center_function > 0d0) & 01110 .or. (lower_function < 0d0 .and. center_function < 0d0)) & 01111 then 01112 lower_volume = center_volume 01113 lower_function = center_function 01114 else 01115 upper_volume = center_volume 01116 upper_function = center_function 01117 end if 01118 end do 01119 01120 new_particle_volume = center_volume 01121 else 01122 ! volume-preserving scheme: Solve the implicit equation: 01123 ! n_part * W(new_vol) * new_vol = weighted_total_volume 01124 ! 01125 ! We assume that the weighting function is strictly monotone 01126 ! so this equation has a unique solution and the solution 01127 ! lies between the min and max particle volumes. We use 01128 ! bisection as this doesn't really need to be fast, just 01129 ! robust. 01130 01131 ! initialize to min and max particle volumes 01132 do i_part = 1,aero_state%bin(i_bin)%n_part 01133 aero_particle => aero_state%bin(i_bin)%particle(i_part) 01134 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01135 aero_data, dry_volume) 01136 if (i_part == 1) then 01137 lower_volume = particle_volume 01138 upper_volume = particle_volume 01139 end if 01140 lower_volume = min(lower_volume, particle_volume) 01141 upper_volume = max(upper_volume, particle_volume) 01142 end do 01143 01144 lower_function = real(aero_state%bin(i_bin)%n_part, kind=dp) 01145 * aero_weight_value(aero_weight, vol2rad(lower_volume)) 01146 * lower_volume - weighted_total_volume 01147 upper_function = real(aero_state%bin(i_bin)%n_part, kind=dp) 01148 * aero_weight_value(aero_weight, vol2rad(upper_volume)) 01149 * upper_volume - weighted_total_volume 01150 01151 ! do 50 rounds of bisection (2^50 = 10^15) 01152 do i_bisect = 1,50 01153 center_volume = (lower_volume + upper_volume) / 2d0 01154 center_function = real(aero_state%bin(i_bin)%n_part, kind=dp) 01155 * aero_weight_value(aero_weight, vol2rad(center_volume)) 01156 * center_volume - weighted_total_volume 01157 if ((lower_function > 0d0 .and. center_function > 0d0) & 01158 .or. (lower_function < 0d0 .and. center_function < 0d0)) & 01159 then 01160 lower_volume = center_volume 01161 lower_function = center_function 01162 else 01163 upper_volume = center_volume 01164 upper_function = center_function 01165 end if 01166 end do 01167 01168 new_particle_volume = center_volume 01169 end if 01170 01171 do i_part = 1,aero_state%bin(i_bin)%n_part 01172 aero_particle => aero_state%bin(i_bin)%particle(i_part) 01173 particle_volume = aero_particle_volume_maybe_dry(aero_particle, & 01174 aero_data, dry_volume) 01175 aero_particle%vol = aero_particle%vol / particle_volume & 01176 * new_particle_volume 01177 if (dry_volume .and. (aero_data%i_water > 0)) then 01178 ! set water to zero if we are doing dry volume averaging 01179 aero_particle%vol(aero_data%i_water) = 0d0 01180 end if 01181 end do 01182 end do 01183 01184 end subroutine aero_state_bin_average_size 01185 01186 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01187 01188 !> Make all particles dry (water set to zero). 01189 subroutine aero_state_make_dry(aero_state, bin_grid, aero_data) 01190 01191 !> Aerosol state to make dry. 01192 type(aero_state_t), intent(inout) :: aero_state 01193 !> Bin grid. 01194 type(bin_grid_t), intent(in) :: bin_grid 01195 !> Aerosol data. 01196 type(aero_data_t), intent(in) :: aero_data 01197 01198 integer :: i_bin, i_part 01199 01200 if (aero_data%i_water > 0) then 01201 do i_bin = 1,bin_grid%n_bin 01202 do i_part = 1,aero_state%bin(i_bin)%n_part 01203 aero_state%bin(i_bin)%particle(i_part)%vol(aero_data%i_water) & 01204 = 0d0 01205 end do 01206 end do 01207 call aero_state_resort(bin_grid, aero_state) 01208 end if 01209 01210 end subroutine aero_state_make_dry 01211 01212 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01213 01214 !> Determines the number of bytes required to pack the given value. 01215 integer function pmc_mpi_pack_size_aero_state(val) 01216 01217 !> Value to pack. 01218 type(aero_state_t), intent(in) :: val 01219 01220 integer :: i, total_size 01221 01222 total_size = 0 01223 total_size = total_size + pmc_mpi_pack_size_real(val%comp_vol) 01224 total_size = total_size + pmc_mpi_pack_size_integer(val%n_part) 01225 total_size = total_size + pmc_mpi_pack_size_integer(size(val%bin)) 01226 do i = 1,size(val%bin) 01227 total_size = total_size + pmc_mpi_pack_size_apa(val%bin(i)) 01228 end do 01229 total_size = total_size + pmc_mpi_pack_size_aia(val%aero_info_array) 01230 pmc_mpi_pack_size_aero_state = total_size 01231 01232 end function pmc_mpi_pack_size_aero_state 01233 01234 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01235 01236 !> Packs the given value into the buffer, advancing position. 01237 subroutine pmc_mpi_pack_aero_state(buffer, position, val) 01238 01239 !> Memory buffer. 01240 character, intent(inout) :: buffer(:) 01241 !> Current buffer position. 01242 integer, intent(inout) :: position 01243 !> Value to pack. 01244 type(aero_state_t), intent(in) :: val 01245 01246 #ifdef PMC_USE_MPI 01247 integer :: prev_position, i 01248 01249 prev_position = position 01250 call pmc_mpi_pack_real(buffer, position, val%comp_vol) 01251 call pmc_mpi_pack_integer(buffer, position, val%n_part) 01252 call pmc_mpi_pack_integer(buffer, position, size(val%bin)) 01253 do i = 1,size(val%bin) 01254 call pmc_mpi_pack_aero_particle_array(buffer, position, val%bin(i)) 01255 end do 01256 call pmc_mpi_pack_aero_info_array(buffer, position, val%aero_info_array) 01257 call assert(850997402, & 01258 position - prev_position <= pmc_mpi_pack_size_aero_state(val)) 01259 #endif 01260 01261 end subroutine pmc_mpi_pack_aero_state 01262 01263 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01264 01265 !> Unpacks the given value from the buffer, advancing position. 01266 subroutine pmc_mpi_unpack_aero_state(buffer, position, val) 01267 01268 !> Memory buffer. 01269 character, intent(inout) :: buffer(:) 01270 !> Current buffer position. 01271 integer, intent(inout) :: position 01272 !> Value to pack. 01273 type(aero_state_t), intent(inout) :: val 01274 01275 #ifdef PMC_USE_MPI 01276 integer :: prev_position, i, n 01277 01278 call aero_state_deallocate(val) 01279 prev_position = position 01280 call pmc_mpi_unpack_real(buffer, position, val%comp_vol) 01281 call pmc_mpi_unpack_integer(buffer, position, val%n_part) 01282 call pmc_mpi_unpack_integer(buffer, position, n) 01283 allocate(val%bin(n)) 01284 do i = 1,size(val%bin) 01285 call aero_particle_array_allocate(val%bin(i)) 01286 call pmc_mpi_unpack_aero_particle_array(buffer, position, val%bin(i)) 01287 end do 01288 call aero_info_array_allocate(val%aero_info_array) 01289 call pmc_mpi_unpack_aero_info_array(buffer, position, val%aero_info_array) 01290 call assert(132104747, & 01291 position - prev_position <= pmc_mpi_pack_size_aero_state(val)) 01292 #endif 01293 01294 end subroutine pmc_mpi_unpack_aero_state 01295 01296 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01297 01298 !> Gathers data from all processes into one aero_state on process 0. 01299 subroutine aero_state_mpi_gather(aero_state, aero_state_total) 01300 01301 !> Local aero_state. 01302 type(aero_state_t), intent(in) :: aero_state 01303 !> Centralized aero_state (only on process 0). 01304 type(aero_state_t), intent(inout) :: aero_state_total 01305 01306 #ifdef PMC_USE_MPI 01307 type(aero_state_t) :: aero_state_transfer 01308 integer :: n_proc, ierr, status(MPI_STATUS_SIZE) 01309 integer :: buffer_size, max_buffer_size, i_proc, position 01310 character, allocatable :: buffer(:) 01311 #endif 01312 01313 if (pmc_mpi_rank() == 0) then 01314 call aero_state_copy(aero_state, aero_state_total) 01315 end if 01316 01317 #ifdef PMC_USE_MPI 01318 01319 if (pmc_mpi_rank() /= 0) then 01320 ! send data from remote processes 01321 max_buffer_size = 0 01322 max_buffer_size = max_buffer_size & 01323 + pmc_mpi_pack_size_aero_state(aero_state) 01324 allocate(buffer(max_buffer_size)) 01325 position = 0 01326 call pmc_mpi_pack_aero_state(buffer, position, aero_state) 01327 call assert(542772170, position <= max_buffer_size) 01328 buffer_size = position 01329 call mpi_send(buffer, buffer_size, MPI_CHARACTER, 0, & 01330 AERO_STATE_TAG_GATHER, MPI_COMM_WORLD, ierr) 01331 call pmc_mpi_check_ierr(ierr) 01332 deallocate(buffer) 01333 else 01334 ! root process receives data from each remote proc 01335 n_proc = pmc_mpi_size() 01336 do i_proc = 1,(n_proc - 1) 01337 ! determine buffer size at root process 01338 call mpi_probe(i_proc, AERO_STATE_TAG_GATHER, MPI_COMM_WORLD, & 01339 status, ierr) 01340 call pmc_mpi_check_ierr(ierr) 01341 call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr) 01342 call pmc_mpi_check_ierr(ierr) 01343 01344 ! get buffer at root process 01345 allocate(buffer(buffer_size)) 01346 call mpi_recv(buffer, buffer_size, MPI_CHARACTER, i_proc, & 01347 AERO_STATE_TAG_GATHER, MPI_COMM_WORLD, status, ierr) 01348 01349 ! unpack it 01350 position = 0 01351 call aero_state_allocate(aero_state_transfer) 01352 call pmc_mpi_unpack_aero_state(buffer, position, & 01353 aero_state_transfer) 01354 call assert(518174881, position == buffer_size) 01355 deallocate(buffer) 01356 01357 call aero_state_add(aero_state_total, aero_state_transfer) 01358 01359 call aero_state_deallocate(aero_state_transfer) 01360 end do 01361 end if 01362 01363 #endif 01364 01365 end subroutine aero_state_mpi_gather 01366 01367 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01368 01369 !> Scatters data from process 0 to all processes by assigning each 01370 !> particle to a random process. 01371 subroutine aero_state_mpi_scatter(aero_state_total, aero_state, & 01372 aero_data) 01373 01374 !> Centralized aero_state (only on process 0. 01375 type(aero_state_t), intent(in) :: aero_state_total 01376 !> Local aero_state. 01377 type(aero_state_t), intent(inout) :: aero_state 01378 !> Aero_data. 01379 type(aero_data_t), intent(in) :: aero_data 01380 01381 #ifndef PMC_USE_MPI 01382 call aero_state_copy(aero_state_total, aero_state) 01383 #else 01384 01385 integer :: n_proc, i_bin, i_part 01386 type(aero_state_t), allocatable :: aero_state_transfers(:) 01387 integer :: ierr, status(MPI_STATUS_SIZE) 01388 integer :: buffer_size, max_buffer_size, i_proc, position 01389 character, allocatable :: buffer(:) 01390 type(aero_particle_t), pointer :: aero_particle 01391 01392 if (pmc_mpi_rank() == 0) then 01393 ! assign particles at random to other processes 01394 allocate(aero_state_transfers(n_proc)) 01395 do i_proc = 1,n_proc 01396 call aero_state_allocate_size(aero_state_transfers(i_proc), & 01397 size(aero_state_total%bin), aero_data%n_spec, aero_data%n_source) 01398 aero_state_transfers(i_proc)%comp_vol = & 01399 aero_state_total%comp_vol / real(n_proc, kind=dp) 01400 end do 01401 do i_bin = 1,size(aero_state_total%bin) 01402 do i_part = 1,aero_state_total%bin(i_bin)%n_part 01403 aero_particle => aero_state_total%bin(i_bin)%particle(i_part) 01404 i_proc = pmc_rand_int(n_proc) 01405 call aero_state_add_particle(aero_state_transfers(i_proc), & 01406 i_bin, aero_particle) 01407 end do 01408 end do 01409 01410 ! just copy our own transfer aero_state 01411 call aero_state_copy(aero_state_transfers(1), aero_state) 01412 01413 ! send the transfer aero_states 01414 n_proc = pmc_mpi_size() 01415 do i_proc = 1,(n_proc - 1) 01416 max_buffer_size = 0 01417 max_buffer_size = max_buffer_size + pmc_mpi_pack_size_aero_state( & 01418 aero_state_transfers(i_proc + 1)) 01419 allocate(buffer(max_buffer_size)) 01420 position = 0 01421 call pmc_mpi_pack_aero_state(buffer, position, & 01422 aero_state_transfers(i_proc + 1)) 01423 call assert(107763122, position <= max_buffer_size) 01424 buffer_size = position 01425 call mpi_send(buffer, buffer_size, MPI_CHARACTER, i_proc, & 01426 AERO_STATE_TAG_SCATTER, MPI_COMM_WORLD, ierr) 01427 call pmc_mpi_check_ierr(ierr) 01428 deallocate(buffer) 01429 end do 01430 01431 ! all done, free data 01432 do i_proc = 1,n_proc 01433 call aero_state_deallocate(aero_state_transfers(i_proc)) 01434 end do 01435 deallocate(aero_state_transfers) 01436 else 01437 ! remote processes just receive data 01438 01439 ! determine buffer size 01440 call mpi_probe(0, AERO_STATE_TAG_SCATTER, MPI_COMM_WORLD, & 01441 status, ierr) 01442 call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr) 01443 call pmc_mpi_check_ierr(ierr) 01444 01445 ! get buffer 01446 allocate(buffer(buffer_size)) 01447 call mpi_recv(buffer, buffer_size, MPI_CHARACTER, 0, & 01448 AERO_STATE_TAG_SCATTER, MPI_COMM_WORLD, status, ierr) 01449 01450 ! unpack it 01451 position = 0 01452 call aero_state_allocate(aero_state) 01453 call pmc_mpi_unpack_aero_state(buffer, position, aero_state) 01454 call assert(374516020, position == buffer_size) 01455 deallocate(buffer) 01456 end if 01457 01458 #endif 01459 01460 end subroutine aero_state_mpi_scatter 01461 01462 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01463 01464 !> Write the aero particle dimension to the given NetCDF file if it 01465 !> is not already present and in any case return the associated 01466 !> dimid. 01467 subroutine aero_state_netcdf_dim_aero_particle(aero_state, ncid, & 01468 dimid_aero_particle) 01469 01470 !> aero_state structure. 01471 type(aero_state_t), intent(in) :: aero_state 01472 !> NetCDF file ID, in data mode. 01473 integer, intent(in) :: ncid 01474 !> Dimid of the aero particle dimension. 01475 integer, intent(out) :: dimid_aero_particle 01476 01477 integer :: status, i_part 01478 integer :: varid_aero_particle 01479 integer :: aero_particle_centers(aero_state%n_part) 01480 01481 ! try to get the dimension ID 01482 status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle) 01483 if (status == NF90_NOERR) return 01484 if (status /= NF90_EBADDIM) call pmc_nc_check(status) 01485 01486 ! dimension not defined, so define now define it 01487 call pmc_nc_check(nf90_redef(ncid)) 01488 01489 call pmc_nc_check(nf90_def_dim(ncid, "aero_particle", & 01490 aero_state%n_part, dimid_aero_particle)) 01491 01492 call pmc_nc_check(nf90_enddef(ncid)) 01493 01494 do i_part = 1,aero_state%n_part 01495 aero_particle_centers(i_part) = i_part 01496 end do 01497 call pmc_nc_write_integer_1d(ncid, aero_particle_centers, & 01498 "aero_particle", (/ dimid_aero_particle /), & 01499 description="dummy dimension variable (no useful value)") 01500 01501 end subroutine aero_state_netcdf_dim_aero_particle 01502 01503 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01504 01505 !> Write the aero removed dimension to the given NetCDF file if it 01506 !> is not already present and in any case return the associated 01507 !> dimid. 01508 subroutine aero_state_netcdf_dim_aero_removed(aero_state, ncid, & 01509 dimid_aero_removed) 01510 01511 !> aero_state structure. 01512 type(aero_state_t), intent(in) :: aero_state 01513 !> NetCDF file ID, in data mode. 01514 integer, intent(in) :: ncid 01515 !> Dimid of the aero removed dimension. 01516 integer, intent(out) :: dimid_aero_removed 01517 01518 integer :: status, i_remove, dim_size 01519 integer :: varid_aero_removed 01520 integer :: aero_removed_centers(max(aero_state%aero_info_array%n_item,1)) 01521 01522 ! try to get the dimension ID 01523 status = nf90_inq_dimid(ncid, "aero_removed", dimid_aero_removed) 01524 if (status == NF90_NOERR) return 01525 if (status /= NF90_EBADDIM) call pmc_nc_check(status) 01526 01527 ! dimension not defined, so define now define it 01528 call pmc_nc_check(nf90_redef(ncid)) 01529 01530 dim_size = max(aero_state%aero_info_array%n_item, 1) 01531 call pmc_nc_check(nf90_def_dim(ncid, "aero_removed", & 01532 dim_size, dimid_aero_removed)) 01533 01534 call pmc_nc_check(nf90_enddef(ncid)) 01535 01536 do i_remove = 1,dim_size 01537 aero_removed_centers(i_remove) = i_remove 01538 end do 01539 call pmc_nc_write_integer_1d(ncid, aero_removed_centers, & 01540 "aero_removed", (/ dimid_aero_removed /), & 01541 description="dummy dimension variable (no useful value)") 01542 01543 end subroutine aero_state_netcdf_dim_aero_removed 01544 01545 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01546 01547 !> Write full state. 01548 subroutine aero_state_output_netcdf(aero_state, ncid, bin_grid, & 01549 aero_data, aero_weight, record_removals, record_optical) 01550 01551 !> aero_state to write. 01552 type(aero_state_t), intent(in) :: aero_state 01553 !> NetCDF file ID, in data mode. 01554 integer, intent(in) :: ncid 01555 !> bin_grid structure. 01556 type(bin_grid_t), intent(in) :: bin_grid 01557 !> aero_data structure. 01558 type(aero_data_t), intent(in) :: aero_data 01559 !> aero_weight structure. 01560 type(aero_weight_t), intent(in) :: aero_weight 01561 !> Whether to output particle removal info. 01562 logical, intent(in) :: record_removals 01563 !> Whether to output aerosol optical properties. 01564 logical, intent(in) :: record_optical 01565 01566 integer :: dimid_aero_particle, dimid_aero_species, dimid_aero_source 01567 integer :: dimid_aero_removed 01568 integer :: i_bin, i_part_in_bin, i_part, i_remove 01569 type(aero_particle_t), pointer :: particle 01570 real(kind=dp) :: aero_particle_mass(aero_state%n_part, aero_data%n_spec) 01571 integer :: aero_n_orig_part(aero_state%n_part, aero_data%n_source) 01572 real(kind=dp) :: aero_absorb_cross_sect(aero_state%n_part) 01573 real(kind=dp) :: aero_scatter_cross_sect(aero_state%n_part) 01574 real(kind=dp) :: aero_asymmetry(aero_state%n_part) 01575 real(kind=dp) :: aero_refract_shell_real(aero_state%n_part) 01576 real(kind=dp) :: aero_refract_shell_imag(aero_state%n_part) 01577 real(kind=dp) :: aero_refract_core_real(aero_state%n_part) 01578 real(kind=dp) :: aero_refract_core_imag(aero_state%n_part) 01579 real(kind=dp) :: aero_core_vol(aero_state%n_part) 01580 integer :: aero_water_hyst_leg(aero_state%n_part) 01581 real(kind=dp) :: aero_comp_vol(aero_state%n_part) 01582 integer :: aero_id(aero_state%n_part) 01583 real(kind=dp) :: aero_least_create_time(aero_state%n_part) 01584 real(kind=dp) :: aero_greatest_create_time(aero_state%n_part) 01585 integer :: aero_removed_id(max(aero_state%aero_info_array%n_item,1)) 01586 integer :: aero_removed_action(max(aero_state%aero_info_array%n_item,1)) 01587 integer :: aero_removed_other_id(max(aero_state%aero_info_array%n_item,1)) 01588 01589 !> \page output_format_aero_state Output File Format: Aerosol Particle State 01590 !! 01591 !! The aerosol state consists of a set of individual aerosol 01592 !! particles, each with its own individual properties. The 01593 !! properties of all particles are stored in arrays, one per 01594 !! property. For example, <tt>aero_absorb_cross_sect(i)</tt> gives 01595 !! the absorption cross section of particle number \c i, while 01596 !! <tt>aero_particle_mass(i,s)</tt> gives the mass of species \c s 01597 !! in particle \c i. The aerosol species are described in \ref 01598 !! output_format_aero_data. 01599 !! 01600 !! Each aerosol particle \c i represents a volume of space known 01601 !! as the computational volume and given by 01602 !! <tt>aero_comp_vol(i)</tt>. Dividing a per-particle quantity by 01603 !! the respective computational volume gives the concentration of 01604 !! that quantity contributed by the particle. For example, summing 01605 !! <tt>aero_particle_mass(i,s) / aero_comp_vol(i)</tt> over all \c 01606 !! i gives the total concentration of species \c s in 01607 !! kg/m^3. Similarly, summing <tt>aero_absorb_cross_sect(i) / 01608 !! aero_comp_vol(i)</tt> over all \c i will give the concentration 01609 !! of scattering cross section in m^2/m^3. 01610 !! 01611 !! The aerosol state uses the \c aero_species NetCDF dimension as 01612 !! specified in the \ref output_format_aero_data section, as well 01613 !! as the NetCDF dimension: 01614 !! - \b aero_particle: number of aerosol particles 01615 !! 01616 !! The aerosol state NetCDF variables are: 01617 !! - \b aero_particle (dim \c aero_particle): dummy dimension variable 01618 !! (no useful value) 01619 !! - \b aero_particle_mass (unit kg, 01620 !! dim <tt>aero_particle x aero_species</tt>): constituent masses of 01621 !! each aerosol particle - <tt>aero_particle_mass(i,s)</tt> gives the 01622 !! mass of species \c s in particle \c i 01623 !! - \b aero_n_orig_part (dim <tt>aero_particle x 01624 !! aero_source</tt>): number of original particles from each 01625 !! source that formed each aerosol particle - 01626 !! <tt>aero_n_orig_part(i,s)</tt> is the number of particles 01627 !! from source \c s that contributed to particle \c i - when 01628 !! particle \c i first enters the simulation (by emissions, 01629 !! dilution, etc.) it has <tt>aero_n_orig_part(i,s) = 1</tt> 01630 !! for the source number \c s it came from (otherwise zero) 01631 !! and when two particles coagulate, their values of \c 01632 !! aero_n_orig_part are added (the number of coagulation 01633 !! events that formed each particle is thus 01634 !! <tt>sum(aero_n_orig_part(i,:)) - 1</tt>) 01635 !! - \b aero_absorb_cross_sect (unit m^2, dim \c aero_particle): 01636 !! optical absorption cross sections of each aerosol particle 01637 !! - \b aero_scatter_cross_sect (unit m^2, dim \c aero_particle): 01638 !! optical scattering cross sections of each aerosol particle 01639 !! - \b aero_asymmetry (dimensionless, dim \c aero_particle): optical 01640 !! asymmetry parameters of each aerosol particle 01641 !! - \b aero_refract_shell_real (dimensionless, dim \c aero_particle): 01642 !! real part of the refractive indices of the shell of each 01643 !! aerosol particle 01644 !! - \b aero_refract_shell_imag (dimensionless, dim \c aero_particle): 01645 !! imaginary part of the refractive indices of the shell of each 01646 !! aerosol particle 01647 !! - \b aero_refract_core_real (dimensionless, dim \c aero_particle): 01648 !! real part of the refractive indices of the core of each 01649 !! aerosol particle 01650 !! - \b aero_refract_core_imag (dimensionless, dim \c aero_particle): 01651 !! imaginary part of the refractive indices of the core of each 01652 !! aerosol particle 01653 !! - \b aero_core_vol (unit m^3, dim \c aero_particle): volume of the 01654 !! optical cores of each aerosol particle 01655 !! - \b aero_water_hyst_leg (dim \c aero_particle): integers 01656 !! specifying which leg of the water hysteresis curve each 01657 !! particle is on, using the MOSAIC numbering convention 01658 !! - \b aero_comp_vol (unit m^3, dim \c aero_particle): computational 01659 !! volume containing each particle 01660 !! - \b aero_id (dim \c aero_particle): unique ID number of each 01661 !! aerosol particle 01662 !! - \b aero_least_create_time (unit s, dim \c aero_particle): least 01663 !! (earliest) creation time of any original constituent particles 01664 !! that coagulated to form each particle, measured from the start 01665 !! of the simulation - a particle is said to be created when it 01666 !! first enters the simulation (by emissions, dilution, etc.) 01667 !! - \b aero_greatest_create_time (unit s, dim \c 01668 !! aero_particle): greatest (latest) creation time of any 01669 !! original constituent particles that coagulated to form each 01670 !! particle, measured from the start of the simulation - a 01671 !! particle is said to be created when it first enters the 01672 !! simulation (by emissions, dilution, etc.) 01673 01674 call aero_data_netcdf_dim_aero_species(aero_data, ncid, & 01675 dimid_aero_species) 01676 call aero_data_netcdf_dim_aero_source(aero_data, ncid, & 01677 dimid_aero_source) 01678 01679 if (aero_state%n_part > 0) then 01680 call aero_state_netcdf_dim_aero_particle(aero_state, ncid, & 01681 dimid_aero_particle) 01682 01683 i_part = 0 01684 do i_bin = 1,bin_grid%n_bin 01685 do i_part_in_bin = 1,aero_state%bin(i_bin)%n_part 01686 i_part = i_part + 1 01687 particle => aero_state%bin(i_bin)%particle(i_part_in_bin) 01688 aero_particle_mass(i_part, :) = particle%vol * aero_data%density 01689 aero_n_orig_part(i_part, :) = particle%n_orig_part 01690 aero_water_hyst_leg(i_part) = particle%water_hyst_leg 01691 aero_comp_vol(i_part) = aero_state%comp_vol & 01692 / aero_weight_value(aero_weight, & 01693 aero_particle_radius(particle)) 01694 aero_id(i_part) = particle%id 01695 aero_least_create_time(i_part) = particle%least_create_time 01696 aero_greatest_create_time(i_part) = particle%greatest_create_time 01697 if (record_optical) then 01698 aero_absorb_cross_sect(i_part) = particle%absorb_cross_sect 01699 aero_scatter_cross_sect(i_part) = particle%scatter_cross_sect 01700 aero_asymmetry(i_part) = particle%asymmetry 01701 aero_refract_shell_real(i_part) = real(particle%refract_shell) 01702 aero_refract_shell_imag(i_part) = & 01703 aimag(particle%refract_shell) 01704 aero_refract_core_real(i_part) = real(particle%refract_core) 01705 aero_refract_core_imag(i_part) = aimag(particle%refract_core) 01706 aero_core_vol(i_part) = particle%core_vol 01707 end if 01708 end do 01709 end do 01710 call pmc_nc_write_real_2d(ncid, aero_particle_mass, & 01711 "aero_particle_mass", (/ dimid_aero_particle, & 01712 dimid_aero_species /), unit="kg", & 01713 long_name="constituent masses of each aerosol particle") 01714 call pmc_nc_write_integer_2d(ncid, aero_n_orig_part, & 01715 "aero_n_orig_part", (/ dimid_aero_particle, & 01716 dimid_aero_source /), & 01717 long_name="number of original constituent particles from " & 01718 // "each source that coagulated to form each aerosol particle") 01719 call pmc_nc_write_integer_1d(ncid, aero_water_hyst_leg, & 01720 "aero_water_hyst_leg", (/ dimid_aero_particle /), & 01721 long_name="leg of the water hysteresis curve leg of each "& 01722 // "aerosol particle") 01723 call pmc_nc_write_real_1d(ncid, aero_comp_vol, & 01724 "aero_comp_vol", (/ dimid_aero_particle /), unit="m^3", & 01725 long_name="computational volume containing each particle") 01726 call pmc_nc_write_integer_1d(ncid, aero_id, & 01727 "aero_id", (/ dimid_aero_particle /), & 01728 long_name="unique ID number of each aerosol particle") 01729 call pmc_nc_write_real_1d(ncid, aero_least_create_time, & 01730 "aero_least_create_time", (/ dimid_aero_particle /), unit="s", & 01731 long_name="least creation time of each aerosol particle", & 01732 description="least (earliest) creation time of any original " & 01733 // "constituent particles that coagulated to form each " & 01734 // "particle, measured from the start of the simulation") 01735 call pmc_nc_write_real_1d(ncid, aero_greatest_create_time, & 01736 "aero_greatest_create_time", (/ dimid_aero_particle /), & 01737 unit="s", & 01738 long_name="greatest creation time of each aerosol particle", & 01739 description="greatest (latest) creation time of any original " & 01740 // "constituent particles that coagulated to form each " & 01741 // "particle, measured from the start of the simulation") 01742 if (record_optical) then 01743 call pmc_nc_write_real_1d(ncid, aero_absorb_cross_sect, & 01744 "aero_absorb_cross_sect", (/ dimid_aero_particle /), & 01745 unit="m^2", & 01746 long_name="optical absorption cross sections of each " & 01747 // "aerosol particle") 01748 call pmc_nc_write_real_1d(ncid, aero_scatter_cross_sect, & 01749 "aero_scatter_cross_sect", (/ dimid_aero_particle /), & 01750 unit="m^2", & 01751 long_name="optical scattering cross sections of each " & 01752 // "aerosol particle") 01753 call pmc_nc_write_real_1d(ncid, aero_asymmetry, & 01754 "aero_asymmetry", (/ dimid_aero_particle /), unit="1", & 01755 long_name="optical asymmetry parameters of each " & 01756 // "aerosol particle") 01757 call pmc_nc_write_real_1d(ncid, aero_refract_shell_real, & 01758 "aero_refract_shell_real", (/ dimid_aero_particle /), & 01759 unit="1", & 01760 long_name="real part of the refractive indices of the " & 01761 // "shell of each aerosol particle") 01762 call pmc_nc_write_real_1d(ncid, aero_refract_shell_imag, & 01763 "aero_refract_shell_imag", (/ dimid_aero_particle /), & 01764 unit="1", & 01765 long_name="imaginary part of the refractive indices of " & 01766 // "the shell of each aerosol particle") 01767 call pmc_nc_write_real_1d(ncid, aero_refract_core_real, & 01768 "aero_refract_core_real", (/ dimid_aero_particle /), & 01769 unit="1", & 01770 long_name="real part of the refractive indices of the core " & 01771 // "of each aerosol particle") 01772 call pmc_nc_write_real_1d(ncid, aero_refract_core_imag, & 01773 "aero_refract_core_imag", (/ dimid_aero_particle /), & 01774 unit="1", & 01775 long_name="imaginary part of the refractive indices of " & 01776 // "the core of each aerosol particle") 01777 call pmc_nc_write_real_1d(ncid, aero_core_vol, & 01778 "aero_core_vol", (/ dimid_aero_particle /), unit="m^3", & 01779 long_name="volume of the optical cores of each " & 01780 // "aerosol particle") 01781 end if 01782 end if 01783 01784 !> \page output_format_aero_removed Output File Format: Aerosol Particle Removal Information 01785 !! 01786 !! When an aerosol particle is introduced into the simulation it 01787 !! is assigned a unique ID number. This ID number will persist 01788 !! over time, allowing tracking of a paticular particle's 01789 !! evolution. If the \c record_removals variable in the input spec 01790 !! file is \c yes, then the every time a particle is removed from 01791 !! the simulation its removal will be recorded in the removal 01792 !! information. 01793 !! 01794 !! The removal information written at timestep \c n contains 01795 !! information about every particle ID that is present at time 01796 !! <tt>(n - 1)</tt> but not present at time \c n. 01797 !! 01798 !! The removal information is always written in the output files, 01799 !! even if no particles were removed in the previous 01800 !! timestep. Unfortunately, NetCDF files cannot contain arrays of 01801 !! length 0. In the case of no particles being removed, the \c 01802 !! aero_removed dimension will be set to 1 and 01803 !! <tt>aero_removed_action(1)</tt> will be 0 (\c AERO_INFO_NONE). 01804 !! 01805 !! When two particles coagulate, the ID number of the combined 01806 !! particle will be the ID particle of the largest constituent, if 01807 !! possible (weighting functions can make this impossible to 01808 !! achieve). A given particle ID may thus be lost due to 01809 !! coagulation (if the resulting combined particle has a different 01810 !! ID), or the ID may be preserved (as the ID of the combined 01811 !! particle). Only if the ID is lost will the particle be recorded 01812 !! in the removal information, and in this case 01813 !! <tt>aero_removed_action(i)</tt> will be 2 (\c AERO_INFO_COAG) 01814 !! and <tt>aero_removed_other_id(i)</tt> will be the ID number of 01815 !! the combined particle. 01816 !! 01817 !! The aerosol removal information NetCDF dimensions are: 01818 !! - \b aero_removed: number of aerosol particles removed from the 01819 !! simulation during the previous timestep (or 1, as described 01820 !! above) 01821 !! 01822 !! The aerosol removal information NetCDF variables are: 01823 !! - \b aero_removed (dim \c aero_removed): dummy dimension variable 01824 !! (no useful value) 01825 !! - \b aero_removed_id (dim \c aero_removed): the ID number of each 01826 !! removed particle 01827 !! - \b aero_removed_action (dim \c aero_removed): the reasons for 01828 !! removal for each particle, with values: 01829 !! - 0 (\c AERO_INFO_NONE): no information (invalid entry) 01830 !! - 1 (\c AERO_INFO_DILUTION): particle was removed due to dilution 01831 !! with outside air 01832 !! - 2 (\c AERO_INFO_COAG): particle was removed due to coagulation 01833 !! - 3 (\c AERO_INFO_HALVED): particle was removed due to halving of 01834 !! the aerosol population 01835 !! - 4 (\c AERO_INFO_WEIGHT): particle was removed due to adjustments 01836 !! in the particle's weighting function 01837 !! - \b aero_removed_other_id (dim \c aero_removed): the ID number of 01838 !! the combined particle formed by coagulation, if the removal reason 01839 !! was coagulation (2, \c AERO_INFO_COAG). May be 0, if the new 01840 !! coagulated particle was not created due to weighting. 01841 01842 if (record_removals) then 01843 call aero_state_netcdf_dim_aero_removed(aero_state, ncid, & 01844 dimid_aero_removed) 01845 if (aero_state%aero_info_array%n_item >= 1) then 01846 do i_remove = 1,aero_state%aero_info_array%n_item 01847 aero_removed_id(i_remove) = & 01848 aero_state%aero_info_array%aero_info(i_remove)%id 01849 aero_removed_action(i_remove) = & 01850 aero_state%aero_info_array%aero_info(i_remove)%action 01851 aero_removed_other_id(i_remove) = & 01852 aero_state%aero_info_array%aero_info(i_remove)%other_id 01853 end do 01854 else 01855 aero_removed_id(1) = 0 01856 aero_removed_action(1) = AERO_INFO_NONE 01857 aero_removed_other_id(1) = 0 01858 end if 01859 call pmc_nc_write_integer_1d(ncid, aero_removed_id, & 01860 "aero_removed_id", (/ dimid_aero_removed /), & 01861 long_name="ID of removed particles") 01862 call pmc_nc_write_integer_1d(ncid, aero_removed_action, & 01863 "aero_removed_action", (/ dimid_aero_removed /), & 01864 long_name="reason for particle removal", & 01865 description="valid is 0 (invalid entry), 1 (removed due to " & 01866 // "dilution), 2 (removed due to coagulation -- combined " & 01867 // "particle ID is in \c aero_removed_other_id), 3 (removed " & 01868 // "due to populating halving), or 4 (removed due to " & 01869 // "weighting changes") 01870 call pmc_nc_write_integer_1d(ncid, aero_removed_other_id, & 01871 "aero_removed_other_id", (/ dimid_aero_removed /), & 01872 long_name="ID of other particle involved in removal", & 01873 description="if <tt>aero_removed_action(i)</tt> is 2 " & 01874 // "(due to coagulation), then " & 01875 // "<tt>aero_removed_other_id(i)</tt> is the ID of the " & 01876 // "resulting combined particle, or 0 if the new particle " & 01877 // "was not created") 01878 end if 01879 01880 end subroutine aero_state_output_netcdf 01881 01882 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 01883 01884 !> Read full state. 01885 subroutine aero_state_input_netcdf(aero_state, ncid, bin_grid, & 01886 aero_data, aero_weight) 01887 01888 !> aero_state to read. 01889 type(aero_state_t), intent(inout) :: aero_state 01890 !> NetCDF file ID, in data mode. 01891 integer, intent(in) :: ncid 01892 !> bin_grid structure. 01893 type(bin_grid_t), intent(in) :: bin_grid 01894 !> aero_data structure. 01895 type(aero_data_t), intent(in) :: aero_data 01896 !> aero_weight structure. 01897 type(aero_weight_t), intent(in) :: aero_weight 01898 01899 integer :: dimid_aero_particle, dimid_aero_removed, n_info_item, n_part 01900 integer :: i_bin, i_part_in_bin, i_part, i_remove, status 01901 type(aero_particle_t) :: aero_particle 01902 character(len=1000) :: name 01903 01904 real(kind=dp), allocatable :: aero_particle_mass(:,:) 01905 integer, allocatable :: aero_n_orig_part(:,:) 01906 real(kind=dp), allocatable :: aero_absorb_cross_sect(:) 01907 real(kind=dp), allocatable :: aero_scatter_cross_sect(:) 01908 real(kind=dp), allocatable :: aero_asymmetry(:) 01909 real(kind=dp), allocatable :: aero_refract_shell_real(:) 01910 real(kind=dp), allocatable :: aero_refract_shell_imag(:) 01911 real(kind=dp), allocatable :: aero_refract_core_real(:) 01912 real(kind=dp), allocatable :: aero_refract_core_imag(:) 01913 real(kind=dp), allocatable :: aero_core_vol(:) 01914 integer, allocatable :: aero_water_hyst_leg(:) 01915 real(kind=dp), allocatable :: aero_comp_vol(:) 01916 integer, allocatable :: aero_id(:) 01917 real(kind=dp), allocatable :: aero_least_create_time(:) 01918 real(kind=dp), allocatable :: aero_greatest_create_time(:) 01919 integer, allocatable :: aero_removed_id(:) 01920 integer, allocatable :: aero_removed_action(:) 01921 integer, allocatable :: aero_removed_other_id(:) 01922 01923 status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle) 01924 if (status == NF90_EBADDIM) then 01925 ! no aero_particle dimension means no particles present 01926 call aero_state_deallocate(aero_state) 01927 call aero_state_allocate_size(aero_state, bin_grid%n_bin, & 01928 aero_data%n_spec, aero_data%n_source) 01929 return 01930 end if 01931 call pmc_nc_check(status) 01932 call pmc_nc_check(nf90_Inquire_Dimension(ncid, dimid_aero_particle, & 01933 name, n_part)) 01934 01935 allocate(aero_particle_mass(n_part, aero_data%n_spec)) 01936 allocate(aero_n_orig_part(n_part, aero_data%n_source)) 01937 allocate(aero_absorb_cross_sect(n_part)) 01938 allocate(aero_scatter_cross_sect(n_part)) 01939 allocate(aero_asymmetry(n_part)) 01940 allocate(aero_refract_shell_real(n_part)) 01941 allocate(aero_refract_shell_imag(n_part)) 01942 allocate(aero_refract_core_real(n_part)) 01943 allocate(aero_refract_core_imag(n_part)) 01944 allocate(aero_core_vol(n_part)) 01945 allocate(aero_water_hyst_leg(n_part)) 01946 allocate(aero_comp_vol(n_part)) 01947 allocate(aero_id(n_part)) 01948 allocate(aero_least_create_time(n_part)) 01949 allocate(aero_greatest_create_time(n_part)) 01950 01951 call pmc_nc_read_real_2d(ncid, aero_particle_mass, & 01952 "aero_particle_mass") 01953 call pmc_nc_read_integer_2d(ncid, aero_n_orig_part, & 01954 "aero_n_orig_part") 01955 call pmc_nc_read_real_1d(ncid, aero_absorb_cross_sect, & 01956 "aero_absorb_cross_sect", must_be_present=.false.) 01957 call pmc_nc_read_real_1d(ncid, aero_scatter_cross_sect, & 01958 "aero_scatter_cross_sect", must_be_present=.false.) 01959 call pmc_nc_read_real_1d(ncid, aero_asymmetry, & 01960 "aero_asymmetry", must_be_present=.false.) 01961 call pmc_nc_read_real_1d(ncid, aero_refract_shell_real, & 01962 "aero_refract_shell_real", must_be_present=.false.) 01963 call pmc_nc_read_real_1d(ncid, aero_refract_shell_imag, & 01964 "aero_refract_shell_imag", must_be_present=.false.) 01965 call pmc_nc_read_real_1d(ncid, aero_refract_core_real, & 01966 "aero_refract_core_real", must_be_present=.false.) 01967 call pmc_nc_read_real_1d(ncid, aero_refract_core_imag, & 01968 "aero_refract_core_imag", must_be_present=.false.) 01969 call pmc_nc_read_real_1d(ncid, aero_core_vol, & 01970 "aero_core_vol", must_be_present=.false.) 01971 call pmc_nc_read_integer_1d(ncid, aero_water_hyst_leg, & 01972 "aero_water_hyst_leg") 01973 call pmc_nc_read_real_1d(ncid, aero_comp_vol, & 01974 "aero_comp_vol") 01975 call pmc_nc_read_integer_1d(ncid, aero_id, & 01976 "aero_id") 01977 call pmc_nc_read_real_1d(ncid, aero_least_create_time, & 01978 "aero_least_create_time") 01979 call pmc_nc_read_real_1d(ncid, aero_greatest_create_time, & 01980 "aero_greatest_create_time") 01981 01982 call aero_state_deallocate(aero_state) 01983 call aero_state_allocate_size(aero_state, bin_grid%n_bin, & 01984 aero_data%n_spec, aero_data%n_source) 01985 01986 call aero_particle_allocate_size(aero_particle, aero_data%n_spec, & 01987 aero_data%n_source) 01988 do i_part = 1,n_part 01989 aero_particle%vol = aero_particle_mass(i_part, :) / aero_data%density 01990 aero_particle%n_orig_part = aero_n_orig_part(i_part, :) 01991 aero_particle%absorb_cross_sect = aero_absorb_cross_sect(i_part) 01992 aero_particle%scatter_cross_sect = aero_scatter_cross_sect(i_part) 01993 aero_particle%asymmetry = aero_asymmetry(i_part) 01994 aero_particle%refract_shell = & 01995 cmplx(aero_refract_shell_real(i_part), & 01996 aero_refract_shell_imag(i_part), kind=dc) 01997 aero_particle%refract_core = cmplx(aero_refract_core_real(i_part), & 01998 aero_refract_core_imag(i_part), kind=dc) 01999 aero_particle%core_vol = aero_core_vol(i_part) 02000 aero_particle%water_hyst_leg = aero_water_hyst_leg(i_part) 02001 aero_state%comp_vol = aero_comp_vol(i_part) & 02002 * aero_weight_value(aero_weight, & 02003 aero_particle_radius(aero_particle)) 02004 aero_particle%id = aero_id(i_part) 02005 aero_particle%least_create_time = aero_least_create_time(i_part) 02006 aero_particle%greatest_create_time = aero_greatest_create_time(i_part) 02007 02008 i_bin = aero_particle_in_bin(aero_particle, bin_grid) 02009 call aero_state_add_particle(aero_state, i_bin, aero_particle) 02010 end do 02011 call aero_particle_deallocate(aero_particle) 02012 02013 deallocate(aero_particle_mass) 02014 deallocate(aero_n_orig_part) 02015 deallocate(aero_absorb_cross_sect) 02016 deallocate(aero_scatter_cross_sect) 02017 deallocate(aero_asymmetry) 02018 deallocate(aero_refract_shell_real) 02019 deallocate(aero_refract_shell_imag) 02020 deallocate(aero_refract_core_real) 02021 deallocate(aero_refract_core_imag) 02022 deallocate(aero_core_vol) 02023 deallocate(aero_water_hyst_leg) 02024 deallocate(aero_comp_vol) 02025 deallocate(aero_id) 02026 deallocate(aero_least_create_time) 02027 deallocate(aero_greatest_create_time) 02028 02029 status = nf90_inq_dimid(ncid, "aero_removed", dimid_aero_removed) 02030 if ((status /= NF90_NOERR) .and. (status /= NF90_EBADDIM)) then 02031 call pmc_nc_check(status) 02032 end if 02033 if (status == NF90_NOERR) then 02034 call pmc_nc_check(nf90_Inquire_Dimension(ncid, dimid_aero_removed, & 02035 name, n_info_item)) 02036 02037 allocate(aero_removed_id(max(n_info_item,1))) 02038 allocate(aero_removed_action(max(n_info_item,1))) 02039 allocate(aero_removed_other_id(max(n_info_item,1))) 02040 02041 call pmc_nc_read_integer_1d(ncid, aero_removed_id, & 02042 "aero_removed_id") 02043 call pmc_nc_read_integer_1d(ncid, aero_removed_action, & 02044 "aero_removed_action") 02045 call pmc_nc_read_integer_1d(ncid, aero_removed_other_id, & 02046 "aero_removed_other_id") 02047 02048 if ((n_info_item > 1) .or. (aero_removed_id(1) /= 0)) then 02049 call aero_info_array_enlarge_to(aero_state%aero_info_array, & 02050 n_info_item) 02051 do i_remove = 1,n_info_item 02052 aero_state%aero_info_array%aero_info(i_remove)%id & 02053 = aero_removed_id(i_remove) 02054 aero_state%aero_info_array%aero_info(i_remove)%action & 02055 = aero_removed_action(i_remove) 02056 aero_state%aero_info_array%aero_info(i_remove)%other_id & 02057 = aero_removed_other_id(i_remove) 02058 end do 02059 end if 02060 02061 deallocate(aero_removed_id) 02062 deallocate(aero_removed_action) 02063 deallocate(aero_removed_other_id) 02064 end if 02065 02066 end subroutine aero_state_input_netcdf 02067 02068 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 02069 02070 end module pmc_aero_state