PartMC 2.1.3
extract_aero_total.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_total program.
00007 
00008 !> Read NetCDF output files and write out the total aerosol number and
00009 !> mass concentrations in text format.
00010 program extract_aero_total
00011 
00012   use netcdf
00013 
00014   integer, parameter :: dp = kind(0.d0)
00015   integer, parameter :: out_unit = 64
00016 
00017   character(len=1000) :: in_prefix, in_filename, out_filename
00018   integer :: ncid
00019   integer :: dimid_aero_species, dimid_aero_particle
00020   integer :: varid_time, varid_aero_species
00021   integer :: varid_aero_particle_mass
00022   integer :: varid_aero_comp_vol
00023   integer :: n_aero_species, n_aero_particle
00024   character(len=1000) :: tmp_str, aero_species_names, remaining_species
00025   real(kind=dp) :: time
00026   real(kind=dp), allocatable :: aero_particle_mass(:,:)
00027   real(kind=dp), allocatable :: aero_comp_vol(:)
00028   integer :: xtype, ndims, nAtts
00029   integer, dimension(nf90_max_var_dims) :: dimids
00030   integer :: ios, i_time, i_part, status
00031   integer :: n_bin, i_bin, n_time, i
00032   real(kind=dp) :: num_conc, mass_conc
00033   character(len=36) :: uuid, run_uuid
00034 
00035   ! process commandline arguments
00036   if (command_argument_count() .ne. 2) then
00037      write(6,*) 'Usage: extract_aero_total <netcdf_state_prefix> ' &
00038           // '<output_filename>'
00039      stop 2
00040   endif
00041   call get_command_argument(1, in_prefix)
00042   call get_command_argument(2, out_filename)
00043 
00044   ! open output file
00045   open(unit=out_unit, file=out_filename, status='replace', iostat=ios)
00046   if (ios /= 0) then
00047      write(0,'(a,a,a,i4)') 'ERROR: unable to open file ', &
00048           trim(out_filename), ' for writing: ', ios
00049      stop 1
00050   end if
00051 
00052   ! write information
00053   write(*,'(a,a)') "Output file: ", trim(out_filename)
00054   write(*,'(a)') "  Each row of output is one time."
00055   write(*,'(a)') "  The columns of output are:"
00056   write(*,'(a)') "    column 1: time (s)"
00057   write(*,'(a)') "    column 2 = aerosol number concentration (#/m^3)"
00058   write(*,'(a)') "    column 3 = aerosol mass concentration (kg/m^3)"
00059 
00060   ! process NetCDF files
00061   i_time = 0
00062   n_time = 0
00063   do while (.true.)
00064      i_time = i_time + 1
00065      write(in_filename,'(a,i8.8,a)') trim(in_prefix), i_time, ".nc"
00066      status = nf90_open(in_filename, NF90_NOWRITE, ncid)
00067      if (status /= NF90_NOERR) then
00068         exit
00069      end if
00070      n_time = i_time
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 concentration = 0
00105         num_conc = 0d0
00106         mass_conc = 0d0
00107      else
00108         call nc_check_msg(status, "getting dimension ID for 'aero_particle'")
00109         call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_aero_particle, &
00110              tmp_str, n_aero_particle), "inquiring dimension 'aero_species'")
00111         
00112         ! read aero_particle_mass
00113         call nc_check_msg(nf90_inq_varid(ncid, "aero_particle_mass", &
00114              varid_aero_particle_mass), &
00115              "getting variable ID for 'aero_particle_mass'")
00116         call nc_check_msg(nf90_Inquire_Variable(ncid, &
00117              varid_aero_particle_mass, tmp_str, xtype, ndims, &
00118              dimids, nAtts), "inquiring variable 'aero_particle_mass'")
00119         if ((ndims /= 2) &
00120              .or. (dimids(1) /= dimid_aero_particle) &
00121              .or. (dimids(2) /= dimid_aero_species)) then
00122            write(0,*) "ERROR: unexpected aero_particle_mass dimids"
00123            stop 1
00124         end if
00125         allocate(aero_particle_mass(n_aero_particle, n_aero_species))
00126         call nc_check_msg(nf90_get_var(ncid, varid_aero_particle_mass, &
00127              aero_particle_mass), "getting variable 'aero_particle_mass'")
00128         
00129         ! read aero_comp_vol
00130         call nc_check_msg(nf90_inq_varid(ncid, "aero_comp_vol", &
00131              varid_aero_comp_vol), "getting variable ID for 'aero_comp_vol'")
00132         call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_comp_vol, &
00133              tmp_str, xtype, ndims, dimids, nAtts), &
00134              "inquiring variable 'aero_comp_vol'")
00135         if ((ndims /= 1) &
00136              .or. (dimids(1) /= dimid_aero_particle)) then
00137            write(0,*) "ERROR: unexpected aero_comp_vol dimids"
00138            stop 1
00139         end if
00140         allocate(aero_comp_vol(n_aero_particle))
00141         call nc_check_msg(nf90_get_var(ncid, varid_aero_comp_vol, &
00142              aero_comp_vol), "getting variable 'aero_comp_vol'")
00143         
00144         call nc_check_msg(nf90_close(ncid), &
00145              "closing file " // trim(in_filename))
00146         
00147         ! compute number and mass concentrations
00148         num_conc = 0d0
00149         mass_conc = 0d0
00150         do i_part = 1,n_aero_particle
00151            num_conc = num_conc + 1d0 / aero_comp_vol(i_part)
00152            mass_conc = mass_conc + sum(aero_particle_mass(i_part,:)) &
00153                 / aero_comp_vol(i_part)
00154         end do
00155         
00156         deallocate(aero_particle_mass)
00157         deallocate(aero_comp_vol)
00158      end if
00159 
00160      ! output data
00161      write(out_unit, '(3e30.15e3)') time, num_conc, mass_conc
00162 
00163   end do
00164 
00165   if (n_time == 0) then
00166      write(*,'(a,a)') 'ERROR: no input file found matching: ', &
00167           trim(in_filename)
00168      stop 1
00169   end if
00170 
00171   close(out_unit)
00172 
00173 contains
00174 
00175 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00176 
00177   !> Check return status of NetCDF function calls.
00178   subroutine nc_check_msg(status, error_msg)
00179 
00180     !> Status return value.
00181     integer, intent(in) :: status
00182     !> Error message in case of failure.
00183     character(len=*), intent(in) :: error_msg
00184 
00185     if (status /= NF90_NOERR) then
00186        write(0,*) trim(error_msg) // " : " // trim(nf90_strerror(status))
00187        stop 1
00188     end if
00189 
00190   end subroutine nc_check_msg
00191 
00192 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00193 
00194 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT
00195   integer function command_argument_count()
00196     command_argument_count = iargc()
00197   end function command_argument_count
00198   subroutine get_command_argument(i, arg)
00199     integer, intent(in) :: i
00200     character(len=*), intent(out) :: arg
00201     call getarg(i, arg)
00202   end subroutine get_command_argument
00203 #endif
00204   
00205 end program extract_aero_total