PartMC 2.1.2
extract_aero_size_mass.F90
Go to the documentation of this file.
00001 ! Copyright (C) 2009-2010 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 extract_aero_size_mass program.
00007 
00008 !> Read NetCDF output files and write out the aerosol mass size
00009 !> distributions in text format.
00010 program extract_aero_size_mass
00011 
00012   use netcdf
00013 
00014   integer, parameter :: dp = kind(0.d0)
00015   integer, parameter :: MAX_N_TIME = 10000
00016   integer, parameter :: out_unit = 64
00017 
00018   character(len=1000) :: in_prefix, in_filename, out_filename
00019   integer :: ncid
00020   integer :: dimid_aero_species, dimid_aero_particle
00021   integer :: varid_time, varid_aero_species
00022   integer :: varid_aero_particle_mass, varid_aero_density
00023   integer :: varid_aero_comp_vol
00024   integer :: n_aero_species, n_aero_particle
00025   character(len=1000) :: tmp_str, aero_species_names
00026   real(kind=dp) :: time
00027   real(kind=dp), allocatable :: aero_particle_mass(:,:)
00028   real(kind=dp), allocatable :: aero_density(:)
00029   real(kind=dp), allocatable :: aero_comp_vol(:)
00030   real(kind=dp), allocatable :: aero_dist(:,:)
00031   integer :: xtype, ndims, nAtts
00032   integer, dimension(nf90_max_var_dims) :: dimids
00033   integer :: ios, i_time, i_spec, i_part, status
00034   integer :: n_bin, i_bin, n_time
00035   real(kind=dp) :: d_min, d_max, diam, volume, log_width
00036   character(len=36) :: uuid, run_uuid
00037 
00038   ! process commandline arguments
00039   if (command_argument_count() .ne. 5) then
00040      write(6,*) 'Usage: extract_aero_size_mass <d_min> <d_max>' &
00041           // ' <n_bin> <netcdf_state_prefix> <output_filename>'
00042      stop 2
00043   endif
00044   call get_command_argument(1, tmp_str)
00045   d_min = string_to_real(tmp_str)
00046   call get_command_argument(2, tmp_str)
00047   d_max = string_to_real(tmp_str)
00048   call get_command_argument(3, tmp_str)
00049   n_bin = string_to_integer(tmp_str)
00050   call get_command_argument(4, in_prefix)
00051   call get_command_argument(5, out_filename)
00052 
00053   ! process NetCDF files
00054   allocate(aero_dist(n_bin, MAX_N_TIME))
00055   aero_dist = 0d0
00056   i_time = 0
00057   n_time = 0
00058   do while (.true.)
00059      i_time = i_time + 1
00060      write(in_filename,'(a,i8.8,a)') trim(in_prefix), i_time, ".nc"
00061      status = nf90_open(in_filename, NF90_NOWRITE, ncid)
00062      if (status /= NF90_NOERR) then
00063         exit
00064      end if
00065      n_time = i_time
00066      if (n_time >= MAX_N_TIME) then
00067         write(0,*) 'ERROR: can only process up to MAX_N_TIME times: ', &
00068              MAX_N_TIME
00069         stop 1
00070      end if
00071 
00072      ! read and check uuid
00073      call nc_check_msg(nf90_get_att(ncid, NF90_GLOBAL, "UUID", uuid), &
00074           "getting global attribute 'UUID'")
00075      if (i_time == 1) then
00076         run_uuid = uuid
00077      else
00078         if (run_uuid /= uuid) then
00079            write(0,*) 'ERROR: UUID mismatch at: ' // trim(in_filename)
00080            stop 1
00081         end if
00082      end if
00083 
00084      ! read time
00085      call nc_check_msg(nf90_inq_varid(ncid, "time", varid_time), &
00086           "getting variable ID for 'time'")
00087      call nc_check_msg(nf90_get_var(ncid, varid_time, time), &
00088           "getting variable 'time'")
00089 
00090      ! read aero_species
00091      call nc_check_msg(nf90_inq_dimid(ncid, "aero_species", &
00092           dimid_aero_species), "getting dimension ID for 'aero_species'")
00093      call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_aero_species, &
00094           tmp_str, n_aero_species), "inquiring dimension 'aero_species'")
00095      call nc_check_msg(nf90_inq_varid(ncid, "aero_species", &
00096           varid_aero_species), "getting variable ID for 'aero_species'")
00097      call nc_check_msg(nf90_get_att(ncid, varid_aero_species, &
00098           "names", aero_species_names), &
00099           "getting attribute 'names' for variable 'aero_species'")
00100      
00101      ! read aero_particle dimension
00102      status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle)
00103      if (status == NF90_EBADDIM) then
00104         ! dimension missing ==> no particles, so skip this time
00105         cycle
00106      end if
00107      call nc_check_msg(status, "getting dimension ID for 'aero_particle'")
00108      call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_aero_particle, &
00109           tmp_str, n_aero_particle), "inquiring dimension 'aero_species'")
00110      
00111      ! read aero_particle_mass
00112      call nc_check_msg(nf90_inq_varid(ncid, "aero_particle_mass", &
00113           varid_aero_particle_mass), &
00114           "getting variable ID for 'aero_particle_mass'")
00115      call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_particle_mass, &
00116           tmp_str, xtype, ndims, dimids, nAtts), &
00117           "inquiring variable 'aero_particle_mass'")
00118      if ((ndims /= 2) &
00119           .or. (dimids(1) /= dimid_aero_particle) &
00120           .or. (dimids(2) /= dimid_aero_species)) then
00121         write(0,*) "ERROR: unexpected aero_particle_mass dimids"
00122         stop 1
00123      end if
00124      allocate(aero_particle_mass(n_aero_particle, n_aero_species))
00125      call nc_check_msg(nf90_get_var(ncid, varid_aero_particle_mass, &
00126           aero_particle_mass), "getting variable 'aero_particle_mass'")
00127      
00128      ! read aero_density
00129      call nc_check_msg(nf90_inq_varid(ncid, "aero_density", &
00130           varid_aero_density), "getting variable ID for 'aero_density'")
00131      call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_density, &
00132           tmp_str, xtype, ndims, dimids, nAtts), &
00133           "inquiring variable 'aero_density'")
00134      if ((ndims /= 1) &
00135           .or. (dimids(1) /= dimid_aero_species)) then
00136         write(0,*) "ERROR: unexpected aero_density dimids"
00137         stop 1
00138      end if
00139      allocate(aero_density(n_aero_species))
00140      call nc_check_msg(nf90_get_var(ncid, varid_aero_density, &
00141           aero_density), "getting variable 'aero_density'")
00142      
00143      ! read aero_comp_vol
00144      call nc_check_msg(nf90_inq_varid(ncid, "aero_comp_vol", &
00145           varid_aero_comp_vol), "getting variable ID for 'aero_comp_vol'")
00146      call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_comp_vol, &
00147           tmp_str, xtype, ndims, dimids, nAtts), &
00148           "inquiring variable 'aero_comp_vol'")
00149      if ((ndims /= 1) &
00150           .or. (dimids(1) /= dimid_aero_particle)) then
00151         write(0,*) "ERROR: unexpected aero_comp_vol dimids"
00152         stop 1
00153      end if
00154      allocate(aero_comp_vol(n_aero_particle))
00155      call nc_check_msg(nf90_get_var(ncid, varid_aero_comp_vol, &
00156           aero_comp_vol), "getting variable 'aero_comp_vol'")
00157      
00158      call nc_check_msg(nf90_close(ncid), &
00159           "closing file " // trim(in_filename))
00160 
00161      ! compute distribution
00162      log_width = (log(d_max) - log(d_min)) / real(n_bin, kind=dp)
00163      do i_part = 1,n_aero_particle
00164         volume = sum(aero_particle_mass(i_part,:) / aero_density)
00165         diam = (volume / (3.14159265358979323846d0 / 6d0))**(1d0/3d0)
00166         i_bin = ceiling((log(diam) - log(d_min)) / log_width)
00167         i_bin = max(1, i_bin)
00168         i_bin = min(n_bin, i_bin)
00169         aero_dist(i_bin, i_time) = aero_dist(i_bin, i_time) &
00170              + sum(aero_particle_mass(i_part,:)) / aero_comp_vol(i_part) &
00171              / log_width
00172      end do
00173 
00174      deallocate(aero_particle_mass)
00175      deallocate(aero_density)
00176      deallocate(aero_comp_vol)
00177   end do
00178 
00179   if (n_time == 0) then
00180      write(*,'(a,a)') 'ERROR: no input file found matching: ', &
00181           trim(in_filename)
00182      stop 1
00183   end if
00184 
00185   ! write information
00186   write(*,'(a,a)') "Output file: ", trim(out_filename)
00187   write(*,'(a)') "  Each row of output is one size bin."
00188   write(*,'(a)') "  The columns of output are:"
00189   write(*,'(a)') "    column   1: bin center diameter (m)"
00190   write(*,'(a)') "    column j+1: mass concentration at time(j) (kg/m^3)"
00191   write(*,'(a)') "  Diameter bins have logarithmic width:"
00192   write(*,'(a,e20.10)') "    log_width = ln(diam(i+1)) - ln(diam(i)) =", &
00193        log_width
00194 
00195   ! open output file
00196   open(unit=out_unit, file=out_filename, status='replace', iostat=ios)
00197   if (ios /= 0) then
00198      write(0,'(a,a,a,i4)') 'ERROR: unable to open file ', &
00199           trim(out_filename), ' for writing: ', ios
00200      stop 1
00201   end if
00202 
00203   ! output data
00204   do i_bin = 1,n_bin
00205      diam = exp((real(i_bin, kind=dp) - 0.5d0) * log_width + log(d_min))
00206      write(out_unit, '(e30.15e3)', advance='no') diam
00207      do i_time = 1,n_time
00208         write(out_unit, '(e30.15e3)', advance='no') &
00209              aero_dist(i_bin, i_time)
00210      end do
00211      write(out_unit, '(a)') ''
00212   end do
00213 
00214   close(out_unit)
00215   deallocate(aero_dist)
00216 
00217 contains
00218 
00219 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00220 
00221   !> Check return status of NetCDF function calls.
00222   subroutine nc_check_msg(status, error_msg)
00223 
00224     !> Status return value.
00225     integer, intent(in) :: status
00226     !> Error message in case of failure.
00227     character(len=*), intent(in) :: error_msg
00228 
00229     if (status /= NF90_NOERR) then
00230        write(0,*) trim(error_msg) // " : " // trim(nf90_strerror(status))
00231        stop 1
00232     end if
00233 
00234   end subroutine nc_check_msg
00235 
00236 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00237 
00238   !> Convert a string to a real.
00239   real(kind=dp) function string_to_real(string)
00240 
00241     !> String to convert.
00242     character(len=*), intent(in) :: string
00243     
00244     real(kind=dp) :: val
00245     integer :: ios
00246 
00247     read(string, '(e40.0)', iostat=ios) val
00248     if (ios /= 0) then
00249        write(0,'(a,a,a,i3)') 'Error converting ', trim(string), &
00250             ' to real: IOSTAT = ', ios
00251        stop 2
00252     end if
00253     string_to_real = val
00254 
00255   end function string_to_real
00256 
00257 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00258 
00259   !> Convert a string to an integer.
00260   integer function string_to_integer(string)
00261 
00262     !> String to convert.
00263     character(len=*), intent(in) :: string
00264     
00265     integer :: val
00266     integer :: ios
00267 
00268     read(string, '(i20)', iostat=ios) val
00269     if (ios /= 0) then
00270        write(0,'(a,a,a,i3)') 'Error converting ', trim(string), &
00271             ' to integer: IOSTAT = ', ios
00272        stop 1
00273     end if
00274     string_to_integer = val
00275 
00276   end function string_to_integer
00277 
00278 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00279 
00280 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT
00281   integer function command_argument_count()
00282     command_argument_count = iargc()
00283   end function command_argument_count
00284   subroutine get_command_argument(i, arg)
00285     integer, intent(in) :: i
00286     character(len=*), intent(out) :: arg
00287     call getarg(i, arg)
00288   end subroutine get_command_argument
00289 #endif
00290   
00291 end program extract_aero_size_mass