PartMC  2.1.5
extract_sectional_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_sectional_aero_size_mass program.
00007 
00008 !> Read NetCDF sectional output files and write out the aerosol mass
00009 !> size distributions in text format.
00010 program extract_sectional_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_diam
00021   integer :: varid_time, varid_aero_species
00022   integer :: varid_aero_diam, varid_aero_diam_widths
00023   integer :: varid_aero_mass_concentration
00024   integer :: n_aero_species, n_bin
00025   character(len=1000) :: tmp_str, aero_species_names
00026   real(kind=dp) :: time
00027   real(kind=dp), allocatable :: aero_dist(:,:)
00028   real(kind=dp), allocatable :: aero_diam(:)
00029   real(kind=dp), allocatable :: aero_diam_widths(:)
00030   real(kind=dp), allocatable :: aero_mass_concentration(:,:)
00031   real(kind=dp), allocatable :: save_aero_diam(:)
00032   integer :: xtype, ndims, nAtts
00033   integer, dimension(nf90_max_var_dims) :: dimids
00034   integer :: ios, i_time, i_spec, i_part, status
00035   integer :: i_bin, n_time, new_n_bin
00036   real(kind=dp) :: log_width
00037   character(len=36) :: uuid, run_uuid
00038 
00039   ! process commandline arguments
00040   if (command_argument_count() .ne. 2) then
00041      write(6,*) 'Usage: extract_sectional_aero_size_mass' &
00042           // ' <netcdf_state_prefix> <output_filename>'
00043      stop 2
00044   endif
00045   call get_command_argument(1, in_prefix)
00046   call get_command_argument(2, out_filename)
00047 
00048   ! process NetCDF files
00049   i_time = 0
00050   n_time = 0
00051   n_bin = 0 ! HACK to shut up gfortran warning
00052   do while (.true.)
00053      i_time = i_time + 1
00054      write(in_filename,'(a,i8.8,a)') trim(in_prefix), i_time, ".nc"
00055      status = nf90_open(in_filename, NF90_NOWRITE, ncid)
00056      if (status /= NF90_NOERR) then
00057         exit
00058      end if
00059      n_time = i_time
00060      if (n_time >= MAX_N_TIME) then
00061         write(0,*) 'ERROR: can only process up to MAX_N_TIME times: ', &
00062              MAX_N_TIME
00063         stop 1
00064      end if
00065 
00066      ! read and check uuid
00067      call nc_check_msg(nf90_get_att(ncid, NF90_GLOBAL, "UUID", uuid), &
00068           "getting global attribute 'UUID'")
00069      if (i_time == 1) then
00070         run_uuid = uuid
00071      else
00072         if (run_uuid /= uuid) then
00073            write(0,*) 'ERROR: UUID mismatch at: ' // trim(in_filename)
00074            stop 1
00075         end if
00076      end if
00077 
00078      ! read time
00079      call nc_check_msg(nf90_inq_varid(ncid, "time", varid_time), &
00080           "getting variable ID for 'time'")
00081      call nc_check_msg(nf90_get_var(ncid, varid_time, time), &
00082           "getting variable 'time'")
00083 
00084      ! read aero_species
00085      call nc_check_msg(nf90_inq_dimid(ncid, "aero_species", &
00086           dimid_aero_species), "getting dimension ID for 'aero_species'")
00087      call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_aero_species, &
00088           tmp_str, n_aero_species), "inquiring dimension 'aero_species'")
00089      call nc_check_msg(nf90_inq_varid(ncid, "aero_species", &
00090           varid_aero_species), "getting variable ID for 'aero_species'")
00091      call nc_check_msg(nf90_get_att(ncid, varid_aero_species, &
00092           "names", aero_species_names), &
00093           "getting attribute 'names' for variable 'aero_species'")
00094      if (i_time == 1) then
00095         write(*,*) "n_aero_species:", n_aero_species
00096         write(*,*) "aero_species_names: ", trim(aero_species_names)
00097      end if
00098      
00099      ! read aero_diam dimension
00100      call nc_check_msg(nf90_inq_dimid(ncid, "aero_diam", &
00101           dimid_aero_diam), "getting dimension ID for 'aero_diam'")
00102      call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_aero_diam, &
00103           tmp_str, new_n_bin), "inquiring dimension 'aero_diam'")
00104 
00105      ! allocate aero_dist on first time
00106      if (i_time == 1) then
00107         n_bin = new_n_bin
00108         allocate(aero_dist(n_bin, MAX_N_TIME))
00109         aero_dist = 0d0
00110      else
00111         if (new_n_bin /= n_bin) then
00112            write(0,*) "ERROR: n_bin changed"
00113            stop 1
00114         end if
00115      end if
00116         
00117      ! read aero_diam variable
00118      call nc_check_msg(nf90_inq_varid(ncid, "aero_diam", &
00119           varid_aero_diam), &
00120           "getting variable ID for 'aero_diam'")
00121      call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_diam, &
00122           tmp_str, xtype, ndims, dimids, nAtts), &
00123           "inquiring variable 'aero_diam'")
00124      if ((ndims /= 1) &
00125           .or. (dimids(1) /= dimid_aero_diam)) then
00126         write(0,*) "ERROR: unexpected aero_diam dimids"
00127         stop 1
00128      end if
00129      allocate(aero_diam(n_bin))
00130      call nc_check_msg(nf90_get_var(ncid, varid_aero_diam, &
00131           aero_diam), "getting variable 'aero_diam'")
00132      if (i_time == 1) then
00133         allocate(save_aero_diam(n_bin))
00134         save_aero_diam = aero_diam
00135      end if
00136 
00137      ! read aero_diam_widths variable
00138      call nc_check_msg(nf90_inq_varid(ncid, "aero_diam_widths", &
00139           varid_aero_diam_widths), &
00140           "getting variable ID 'aero_diam_widths'")
00141      call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_diam_widths, &
00142           tmp_str, xtype, ndims, dimids, nAtts), &
00143           "inquiring variable 'aero_diam_widths'")
00144      if ((ndims /= 1) &
00145           .or. (dimids(1) /= dimid_aero_diam)) then
00146         write(0,*) "ERROR: unexpected aero_diam_widths dimids"
00147         stop 1
00148      end if
00149      allocate(aero_diam_widths(n_bin))
00150      call nc_check_msg(nf90_get_var(ncid, varid_aero_diam_widths, &
00151           aero_diam_widths), "getting variable 'aero_diam_widths'")
00152 
00153      ! read aero_mass_concentration
00154      call nc_check_msg(nf90_inq_varid(ncid, "aero_mass_concentration", &
00155           varid_aero_mass_concentration), &
00156           "getting variable ID for 'aero_mass_concentration'")
00157      call nc_check_msg(nf90_Inquire_Variable(ncid, &
00158           varid_aero_mass_concentration, tmp_str, xtype, ndims, &
00159           dimids, nAtts), "inquiring variable 'aero_mass_concentration'")
00160      if ((ndims /= 2) &
00161           .or. (dimids(1) /= dimid_aero_diam) &
00162           .or. (dimids(2) /= dimid_aero_species)) then
00163         write(0,*) "ERROR: unexpected aero_mass_concentration dimids"
00164         stop 1
00165      end if
00166      allocate(aero_mass_concentration(n_bin, n_aero_species))
00167      call nc_check_msg(nf90_get_var(ncid, varid_aero_mass_concentration, &
00168           aero_mass_concentration), &
00169           "getting variable 'aero_mass_concentration'")
00170      
00171      call nc_check_msg(nf90_close(ncid), &
00172           "closing file " // trim(in_filename))
00173 
00174      ! compute distribution
00175      log_width = aero_diam_widths(1)
00176      do i_bin = 1,n_bin
00177         aero_dist(i_bin, i_time) = aero_dist(i_bin, i_time) &
00178              + sum(aero_mass_concentration(i_bin,:))
00179      end do
00180 
00181      deallocate(aero_diam)
00182      deallocate(aero_diam_widths)
00183      deallocate(aero_mass_concentration)
00184   end do
00185 
00186   if (n_time == 0) then
00187      write(*,'(a,a)') 'ERROR: no input file found matching: ', &
00188           trim(in_filename)
00189      stop 1
00190   end if
00191 
00192   ! write information
00193   write(*,*) "Output file array A has:"
00194   write(*,*) "  A(i, 1) = diameter(i) (m)"
00195   write(*,*) "  A(i, j+1) = mass concentration at diameter(i) and " &
00196        // "time(j) (kg/m^3)"
00197   write(*,*) "Diameter bins have logarithmic width:"
00198   write(*,*) "  log_width = ln(diameter(i+1)) - ln(diameter(i)) =", log_width
00199 
00200   ! open output file
00201   open(unit=out_unit, file=out_filename, status='replace', iostat=ios)
00202   if (ios /= 0) then
00203      write(0,'(a,a,a,i4)') 'ERROR: unable to open file ', &
00204           trim(out_filename), ' for writing: ', ios
00205      stop 1
00206   end if
00207 
00208   ! output data
00209   do i_bin = 1,n_bin
00210      write(out_unit, '(e30.15e3)', advance='no') &
00211           save_aero_diam(i_bin)
00212      do i_time = 1,n_time
00213         write(out_unit, '(e30.15e3)', advance='no') &
00214              aero_dist(i_bin, i_time)
00215      end do
00216      write(out_unit, '(a)') ''
00217   end do
00218 
00219   close(out_unit)
00220   deallocate(save_aero_diam)
00221   deallocate(aero_dist)
00222 
00223 contains
00224 
00225 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00226 
00227   !> Check return status of NetCDF function calls.
00228   subroutine nc_check_msg(status, error_msg)
00229 
00230     !> Status return value.
00231     integer, intent(in) :: status
00232     !> Error message in case of failure.
00233     character(len=*), intent(in) :: error_msg
00234 
00235     if (status /= NF90_NOERR) then
00236        write(0,*) trim(error_msg) // " : " // trim(nf90_strerror(status))
00237        stop 1
00238     end if
00239 
00240   end subroutine nc_check_msg
00241 
00242 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00243 
00244 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT
00245   integer function command_argument_count()
00246     command_argument_count = iargc()
00247   end function command_argument_count
00248   subroutine get_command_argument(i, arg)
00249     integer, intent(in) :: i
00250     character(len=*), intent(out) :: arg
00251     call getarg(i, arg)
00252   end subroutine get_command_argument
00253 #endif
00254   
00255 end program extract_sectional_aero_size_mass