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