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