PartMC 2.1.3
bin_grid.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_bin_grid module.
00007 
00008 !> The bin_grid_t structure and associated subroutines.
00009 module pmc_bin_grid
00010 
00011   use pmc_constants
00012   use pmc_util
00013   use pmc_spec_file
00014   use pmc_mpi
00015   use pmc_netcdf
00016 #ifdef PMC_USE_MPI
00017   use mpi
00018 #endif
00019 
00020   !> 1D grid of size bins.
00021   !!
00022   !! The grid of bins is logarithmically spaced in volume, an
00023   !! assumption that is quite heavily incorporated into the code. At
00024   !! some point in the future it would be nice to relax this
00025   !! assumption.
00026   type bin_grid_t
00027      !> Number of bins.
00028      integer :: n_bin
00029      !> Len n_bin, bin center radii (m^3).
00030      real(kind=dp), pointer :: center_radius(:)
00031      !> Len (n_bin + 1), bin edge radii (m^3).
00032      real(kind=dp), pointer :: edge_radius(:)
00033      !> Bin logarithmic width, equal to <tt>log(edge_radius(i+1)) -
00034      !> log(edge_radius(i))</tt> for any \c i (dimensionless).
00035      real(kind=dp) :: log_width
00036   end type bin_grid_t
00037 
00038 contains
00039   
00040 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00041 
00042   !> Allocates a bin_grid.
00043   subroutine bin_grid_allocate(bin_grid)
00044 
00045     !> Bin grid.
00046     type(bin_grid_t), intent(out) :: bin_grid
00047 
00048     bin_grid%n_bin = 0
00049     allocate(bin_grid%center_radius(0))
00050     allocate(bin_grid%edge_radius(0))
00051 
00052   end subroutine bin_grid_allocate
00053 
00054 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00055 
00056   !> Allocates a bin_grid of the given size.
00057   subroutine bin_grid_allocate_size(bin_grid, n_bin)
00058 
00059     !> Bin grid.
00060     type(bin_grid_t), intent(out) :: bin_grid
00061     !> Number of bins.
00062     integer, intent(in) :: n_bin
00063 
00064     bin_grid%n_bin = n_bin
00065     allocate(bin_grid%center_radius(n_bin))
00066     allocate(bin_grid%edge_radius(n_bin + 1))
00067 
00068   end subroutine bin_grid_allocate_size
00069 
00070 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00071 
00072   !> Frees all memory.
00073   subroutine bin_grid_deallocate(bin_grid)
00074 
00075     !> Bin_grid to free.
00076     type(bin_grid_t), intent(inout) :: bin_grid
00077 
00078     deallocate(bin_grid%center_radius)
00079     deallocate(bin_grid%edge_radius)
00080 
00081   end subroutine bin_grid_deallocate
00082 
00083 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00084 
00085   !> Convert a concentration f(vol)d(vol) to f(ln(r))d(ln(r))
00086   !> where vol = 4/3 pi r^3.
00087   subroutine vol_to_lnr(r, f_vol, f_lnr)
00088     
00089     !> Radius (m).
00090     real(kind=dp), intent(in) :: r
00091     !> Concentration as a function of volume.
00092     real(kind=dp), intent(in) :: f_vol
00093     !> Concentration as a function of ln(r).
00094     real(kind=dp), intent(out) :: f_lnr
00095     
00096     f_lnr = f_vol * 4d0 * const%pi * r**3
00097     
00098   end subroutine vol_to_lnr
00099   
00100 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00101 
00102   !> Generates the bin grid given the range and number of bins.
00103   subroutine bin_grid_make(bin_grid, n_bin, r_min, r_max)
00104     
00105     !> New bin grid.
00106     type(bin_grid_t), intent(inout) :: bin_grid
00107     !> Number of bins.
00108     integer, intent(in) :: n_bin
00109     !> Minimum radius (m^3).
00110     real(kind=dp), intent(in) :: r_min
00111     !> Minimum radius (m^3).
00112     real(kind=dp), intent(in) :: r_max
00113 
00114     integer :: i_bin
00115 
00116     call assert_msg(538534122, n_bin > 0, &
00117          "bin_grid requires a positive n_bin, not: " &
00118          // trim(integer_to_string(n_bin)))
00119     call assert_msg(966541762, r_min > 0d0, &
00120          "bin_grid requires a positive r_min, not: " &
00121          // trim(real_to_string(r_min)))
00122     call assert_msg(966541762, r_min < r_max, &
00123          "bin_grid requires r_min < r_max, not: " &
00124          // trim(real_to_string(r_min)) // " and " &
00125          // trim(real_to_string(r_max)))
00126     call bin_grid_deallocate(bin_grid)
00127     call bin_grid_allocate_size(bin_grid, n_bin)
00128     call logspace(r_min, r_max, bin_grid%edge_radius)
00129     do i_bin = 1,n_bin
00130        bin_grid%center_radius(i_bin) &
00131             = exp(0.5d0 * log(bin_grid%edge_radius(i_bin)) &
00132             + 0.5d0 * log(bin_grid%edge_radius(i_bin + 1)))
00133     end do
00134     bin_grid%log_width = (log(r_max) - log(r_min)) / real(n_bin, kind=dp)
00135 
00136   end subroutine bin_grid_make
00137 
00138 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00139 
00140   !> Find the bin number that contains a given particle.
00141   !!
00142   !! This assumes logarithmically spaced bins. If a particle is below
00143   !! the smallest bin or above the largest bin, then it is returned as
00144   !! being in the smallest or largest bin, respectively.
00145   integer function bin_grid_particle_in_bin(bin_grid, radius)
00146 
00147     !> Bin_grid.
00148     type(bin_grid_t), intent(in) :: bin_grid
00149     !> Radius of particle.
00150     real(kind=dp), intent(in) :: radius
00151 
00152     call assert(448215689, bin_grid%n_bin >= 1)
00153     bin_grid_particle_in_bin = logspace_find(bin_grid%edge_radius(1), &
00154          bin_grid%edge_radius(bin_grid%n_bin + 1), bin_grid%n_bin + 1, &
00155          radius)
00156 
00157   end function bin_grid_particle_in_bin
00158   
00159 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00160 
00161   !> Read the specification for a bin_grid from a spec file and
00162   !> generate it.
00163   subroutine spec_file_read_bin_grid(file, bin_grid)
00164 
00165     !> Spec file.
00166     type(spec_file_t), intent(inout) :: file
00167     !> Bin grid.
00168     type(bin_grid_t), intent(inout) :: bin_grid
00169 
00170     integer :: n_bin
00171     real(kind=dp) :: d_min, d_max
00172 
00173     !> \page input_format_bin_grid Input File Format: Diameter Axis Bin Grid
00174     !!
00175     !! The bin grid is logarithmic in diameter, consisting of
00176     !! \f$n_{\rm bin}\f$ bins with centers \f$c_i\f$ (\f$i =
00177     !! 1,\ldots,n_{\rm bin}\f$) and edges \f$e_i\f$ (\f$i =
00178     !! 1,\ldots,(n_{\rm bin} + 1)\f$) such that \f$e_{i+1}/e_i\f$ is a
00179     !! constant and \f$c_i/e_i = \sqrt{e_{i+1}/e_i}\f$. That is,
00180     !! \f$\ln(e_i)\f$ are uniformly spaced and \f$\ln(c_i)\f$ are the
00181     !! arithmetic centers.
00182     !!
00183     !! The diameter axis bin grid is specified by the parameters:
00184     !!   - \b n_bin (integer): The number of bins \f$n_{\rm bin}\f$.
00185     !!   - \b d_min (real, unit m): The left edge of the left-most bin,
00186     !!     \f$e_1\f$.
00187     !!   - \b d_max (real, unit m): The right edge of the right-most bin,
00188     !!     \f$e_{n_{\rm bin} + 1}\f$.
00189     !!
00190     !! See also:
00191     !!   - \ref spec_file_format --- the input file text format
00192     !!   - \ref output_format_bin_grid --- the corresponding output format
00193 
00194     call spec_file_read_integer(file, 'n_bin', n_bin)
00195     call spec_file_read_real(file, 'd_min', d_min)
00196     call spec_file_read_real(file, 'd_max', d_max)
00197     call bin_grid_make(bin_grid, n_bin, diam2rad(d_min), diam2rad(d_max))
00198 
00199   end subroutine spec_file_read_bin_grid
00200 
00201 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00202 
00203   !> Determines the number of bytes required to pack the given value.
00204   integer function pmc_mpi_pack_size_bin_grid(val)
00205 
00206     !> Value to pack.
00207     type(bin_grid_t), intent(in) :: val
00208 
00209     pmc_mpi_pack_size_bin_grid = &
00210          pmc_mpi_pack_size_integer(val%n_bin) &
00211          + pmc_mpi_pack_size_real_array(val%center_radius) &
00212          + pmc_mpi_pack_size_real_array(val%edge_radius) &
00213          + pmc_mpi_pack_size_real(val%log_width)
00214 
00215   end function pmc_mpi_pack_size_bin_grid
00216 
00217 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00218 
00219   !> Packs the given value into the buffer, advancing position.
00220   subroutine pmc_mpi_pack_bin_grid(buffer, position, val)
00221 
00222     !> Memory buffer.
00223     character, intent(inout) :: buffer(:)
00224     !> Current buffer position.
00225     integer, intent(inout) :: position
00226     !> Value to pack.
00227     type(bin_grid_t), intent(in) :: val
00228 
00229 #ifdef PMC_USE_MPI
00230     integer :: prev_position
00231 
00232     prev_position = position
00233     call pmc_mpi_pack_integer(buffer, position, val%n_bin)
00234     call pmc_mpi_pack_real_array(buffer, position, val%center_radius)
00235     call pmc_mpi_pack_real_array(buffer, position, val%edge_radius)
00236     call pmc_mpi_pack_real(buffer, position, val%log_width)
00237     call assert(385455586, &
00238          position - prev_position <= pmc_mpi_pack_size_bin_grid(val))
00239 #endif
00240 
00241   end subroutine pmc_mpi_pack_bin_grid
00242 
00243 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00244 
00245   !> Unpacks the given value from the buffer, advancing position.
00246   subroutine pmc_mpi_unpack_bin_grid(buffer, position, val)
00247 
00248     !> Memory buffer.
00249     character, intent(inout) :: buffer(:)
00250     !> Current buffer position.
00251     integer, intent(inout) :: position
00252     !> Value to pack.
00253     type(bin_grid_t), intent(inout) :: val
00254 
00255 #ifdef PMC_USE_MPI
00256     integer :: prev_position
00257 
00258     prev_position = position
00259     call pmc_mpi_unpack_integer(buffer, position, val%n_bin)
00260     call pmc_mpi_unpack_real_array(buffer, position, val%center_radius)
00261     call pmc_mpi_unpack_real_array(buffer, position, val%edge_radius)
00262     call pmc_mpi_unpack_real(buffer, position, val%log_width)
00263     call assert(741838730, &
00264          position - prev_position <= pmc_mpi_pack_size_bin_grid(val))
00265 #endif
00266 
00267   end subroutine pmc_mpi_unpack_bin_grid
00268 
00269 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00270 
00271   !> Write the aero_diam dimension to the given NetCDF file if it is
00272   !> not already present and in any case return the associated dimid.
00273   subroutine bin_grid_netcdf_dim_aero_diam(bin_grid, ncid, &
00274        dimid_aero_diam)
00275 
00276     !> Bin_grid structure.
00277     type(bin_grid_t), intent(in) :: bin_grid
00278     !> NetCDF file ID, in data mode.
00279     integer, intent(in) :: ncid
00280     !> Dimid of the aero_diam dimension.
00281     integer, intent(out) :: dimid_aero_diam
00282 
00283     integer :: status, varid_aero_diam
00284     integer :: dimid_aero_diam_edges, varid_aero_diam_edges, 
00285          varid_aero_diam_widths
00286     real(kind=dp) :: aero_diam_centers(bin_grid%n_bin), 
00287          aero_diam_edges(bin_grid%n_bin + 1)
00288     real(kind=dp) :: aero_diam_widths(bin_grid%n_bin)
00289 
00290     status = nf90_inq_dimid(ncid, "aero_diam", dimid_aero_diam)
00291     if (status == NF90_NOERR) return
00292     if (status /= NF90_EBADDIM) call pmc_nc_check(status)
00293 
00294     ! dimension not defined, so define now define it
00295     call pmc_nc_check(nf90_redef(ncid))
00296 
00297     call pmc_nc_check(nf90_def_dim(ncid, "aero_diam", &
00298          bin_grid%n_bin, dimid_aero_diam))
00299     call pmc_nc_check(nf90_def_dim(ncid, "aero_diam_edges", &
00300          bin_grid%n_bin + 1, dimid_aero_diam_edges))
00301 
00302     call pmc_nc_check(nf90_enddef(ncid))
00303 
00304     aero_diam_centers = rad2diam(bin_grid%center_radius)
00305     aero_diam_widths = bin_grid%log_width
00306     aero_diam_edges = rad2diam(bin_grid%edge_radius)
00307 
00308     call pmc_nc_write_real_1d(ncid, aero_diam_centers, &
00309          "aero_diam", (/ dimid_aero_diam /), unit="m", &
00310          long_name="aerosol radius axis bin centers", &
00311          description="logarithmically spaced centers of radius axis grid, " &
00312          // "so that aero_diam(i) / aero_diam_edges(i) = " &
00313          // "0.5 * aero_diam_edges(i+1) / aero_diam_edges(i)")
00314     call pmc_nc_write_real_1d(ncid, aero_diam_edges, &
00315          "aero_diam_edges", (/ dimid_aero_diam_edges /), unit="m", &
00316          long_name="aerosol radius axis bin edges", &
00317          description="logarithmically spaced edges of radius axis grid, " &
00318          // "with one more edge than center")
00319     call pmc_nc_write_real_1d(ncid, aero_diam_widths, &
00320          "aero_diam_widths", (/ dimid_aero_diam /), unit="m", &
00321          long_name="aerosol radius axis bin widths", &
00322          description="base-e logarithmic widths of radius axis grid, " &
00323          // "so that aero_diam_widths(i)" &
00324          // "= ln(aero_diam_edges(i+1) / aero_diam_edges(i)) and " &
00325          // "all bins have the same width")
00326 
00327   end subroutine bin_grid_netcdf_dim_aero_diam
00328 
00329 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00330 
00331   !> Write full state.
00332   subroutine bin_grid_output_netcdf(bin_grid, ncid)
00333     
00334     !> bin_grid to write.
00335     type(bin_grid_t), intent(in) :: bin_grid
00336     !> NetCDF file ID, in data mode.
00337     integer, intent(in) :: ncid
00338 
00339     integer :: dimid_aero_diam
00340 
00341     !> \page output_format_bin_grid Output File Format: Bin Grid Data
00342     !!
00343     !! The bin grid data NetCDF dimensions are:
00344     !!   - \b aero_diam: number of bins (grid cells) on the diameter axis
00345     !!   - \b aero_diam_edges: number of bin edges (grid cell edges) on
00346     !!     the diameter axis --- always equal to <tt>aero_diam + 1</tt>
00347     !!
00348     !! The bin grid data NetCDF variables are:
00349     !!   - \b aero_diam (unit m, dim \c aero_diam): aerosol diameter axis
00350     !!     bin centers --- centered on a logarithmic scale from the edges, so
00351     !!     that <tt>aero_diam(i) / aero_diam_edges(i) =
00352     !!     sqrt(aero_diam_edges(i+1) / aero_diam_edges(i))</tt>
00353     !!   - \b aero_diam_edges (unit m, dim \c aero_diam_edges): aersol
00354     !!     diameter axis bin edges (there is one more edge than center)
00355     !!   - \b aero_diam_widths (dimensionless, dim \c aero_diam):
00356     !!     the base-e logarithmic bin widths --- <tt>aero_diam_widths(i)
00357     !!     = ln(aero_diam_edges(i+1) / aero_diam_edges(i))</tt>, so
00358     !!     all bins have the same width
00359     !!
00360     !! See also:
00361     !!   - \ref input_format_bin_grid --- the corresponding input format
00362 
00363     call bin_grid_netcdf_dim_aero_diam(bin_grid, ncid, &
00364          dimid_aero_diam)
00365 
00366     ! no need to write any more data as it's all contained in the
00367     ! dimension and associated variables
00368 
00369   end subroutine bin_grid_output_netcdf
00370 
00371 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00372 
00373   !> Read full state.
00374   subroutine bin_grid_input_netcdf(bin_grid, ncid)
00375     
00376     !> bin_grid to read.
00377     type(bin_grid_t), intent(inout) :: bin_grid
00378     !> NetCDF file ID, in data mode.
00379     integer, intent(in) :: ncid
00380 
00381     integer :: dimid_aero_diam, n_bin
00382     character(len=1000) :: name
00383     real(kind=dp), allocatable :: aero_diam_centers(:)
00384     real(kind=dp), allocatable :: aero_diam_edges(:)
00385     real(kind=dp), allocatable :: aero_diam_widths(:)
00386 
00387     call pmc_nc_check(nf90_inq_dimid(ncid, "aero_diam", dimid_aero_diam))
00388     call pmc_nc_check(nf90_Inquire_Dimension(ncid, dimid_aero_diam, name, &
00389          n_bin))
00390 
00391     call bin_grid_deallocate(bin_grid)
00392     call bin_grid_allocate_size(bin_grid, n_bin)
00393 
00394     allocate(aero_diam_centers(n_bin))
00395     allocate(aero_diam_edges(n_bin + 1))
00396     allocate(aero_diam_widths(n_bin))
00397 
00398     call pmc_nc_read_real_1d(ncid, aero_diam_centers, &
00399          "aero_diam")
00400     call pmc_nc_read_real_1d(ncid, aero_diam_edges, &
00401          "aero_diam_edges")
00402     call pmc_nc_read_real_1d(ncid, aero_diam_widths, &
00403          "aero_diam_widths")
00404 
00405     bin_grid%center_radius = diam2rad(aero_diam_centers)
00406     bin_grid%edge_radius = diam2rad(aero_diam_edges)
00407     bin_grid%log_width = aero_diam_widths(1)
00408 
00409     deallocate(aero_diam_centers)
00410     deallocate(aero_diam_edges)
00411     deallocate(aero_diam_widths)
00412 
00413   end subroutine bin_grid_input_netcdf
00414 
00415 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00416   
00417 end module pmc_bin_grid