PartMC 2.1.4
extract_aero_species.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_species program.
00007 
00008 !> Read NetCDF output files and write out the per-species bulk aerosol
00009 !> concentrations in text format.
00010 program extract_aero_species
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   real(kind=dp), allocatable :: aero_conc(:)
00029   integer :: xtype, ndims, nAtts
00030   integer, dimension(nf90_max_var_dims) :: dimids
00031   integer :: ios, i_time, i_spec, i_part, status, n_time, i
00032   character(len=36) :: uuid, run_uuid
00033 
00034   ! process commandline arguments
00035   if (command_argument_count() .ne. 2) then
00036      write(6,*) 'Usage: extract_aero_species <netcdf_state_prefix>' &
00037           // ' <output_filename>'
00038      stop 2
00039   endif
00040   call get_command_argument(1, in_prefix)
00041   call get_command_argument(2, out_filename)
00042 
00043   ! open output file
00044   open(unit=out_unit, file=out_filename, status='replace', iostat=ios)
00045   if (ios /= 0) then
00046      write(0,'(a,a,a,i4)') 'ERROR: unable to open file ', &
00047           trim(out_filename), ' for writing: ', ios
00048      stop 1
00049   end if
00050 
00051   ! write information
00052   write(*,'(a,a)') "Output file: ", trim(out_filename)
00053   write(*,'(a)') "  Each row of output is one time."
00054   write(*,'(a)') "  The columns of output are:"
00055   write(*,'(a)') "    column  1: time (s)"
00056 
00057   ! process NetCDF files
00058   i_time = 0
00059   n_time = 0
00060   do while (.true.)
00061      i_time = i_time + 1
00062      write(in_filename,'(a,i8.8,a)') trim(in_prefix), i_time, ".nc"
00063      status = nf90_open(in_filename, NF90_NOWRITE, ncid)
00064      if (status /= NF90_NOERR) then
00065         exit
00066      end if
00067      n_time = i_time
00068 
00069      ! read and check uuid
00070      call nc_check_msg(nf90_get_att(ncid, NF90_GLOBAL, "UUID", uuid), &
00071           "getting global attribute 'UUID'")
00072      if (i_time == 1) then
00073         run_uuid = uuid
00074      else
00075         if (run_uuid /= uuid) then
00076            write(0,*) 'ERROR: UUID mismatch at: ' // trim(in_filename)
00077            stop 1
00078         end if
00079      end if
00080 
00081      ! read time
00082      call nc_check_msg(nf90_inq_varid(ncid, "time", varid_time), &
00083           "getting variable ID for 'time'")
00084      call nc_check_msg(nf90_get_var(ncid, varid_time, time), &
00085           "getting variable 'time'")
00086 
00087      ! read aero_species
00088      call nc_check_msg(nf90_inq_dimid(ncid, "aero_species", &
00089           dimid_aero_species), "getting dimension ID for 'aero_species'")
00090      call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_aero_species, &
00091           tmp_str, n_aero_species), "inquiring dimension 'aero_species'")
00092      call nc_check_msg(nf90_inq_varid(ncid, "aero_species", &
00093           varid_aero_species), "getting variable ID for 'aero_species'")
00094      call nc_check_msg(nf90_get_att(ncid, varid_aero_species, &
00095           "names", aero_species_names), &
00096           "getting attribute 'names' for variable 'aero_species'")
00097      if (i_time == 1) then
00098         allocate(aero_conc(n_aero_species))
00099         remaining_species = aero_species_names
00100         do i_spec = 1,n_aero_species
00101            if (i_spec < n_aero_species) then
00102               i = index(remaining_species, ',')
00103            else
00104               i = index(remaining_species, ' ')
00105            end if
00106            write(*,'(a,i2,a,a,a)') '    column ', i_spec + 1, ': ', &
00107                 remaining_species(1:(i-1)), ' (kg/m^3)'
00108            remaining_species = remaining_species((i+1):)
00109         end do
00110      end if
00111      
00112      ! read aero_particle dimension
00113      status = nf90_inq_dimid(ncid, "aero_particle", dimid_aero_particle)
00114      if (status == NF90_EBADDIM) then
00115         ! dimension missing ==> no particles, so concentration = 0
00116         aero_conc = 0d0
00117      else
00118         call nc_check_msg(status, "getting dimension ID for 'aero_particle'")
00119         call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_aero_particle, &
00120              tmp_str, n_aero_particle), "inquiring dimension 'aero_species'")
00121         
00122         ! read aero_particle_mass
00123         call nc_check_msg(nf90_inq_varid(ncid, "aero_particle_mass", &
00124              varid_aero_particle_mass), &
00125              "getting variable ID for 'aero_particle_mass'")
00126         call nc_check_msg(nf90_Inquire_Variable(ncid, &
00127              varid_aero_particle_mass, tmp_str, xtype, ndims, &
00128              dimids, nAtts), "inquiring variable 'aero_particle_mass'")
00129         if ((ndims /= 2) &
00130              .or. (dimids(1) /= dimid_aero_particle) &
00131              .or. (dimids(2) /= dimid_aero_species)) then
00132            write(0,*) "ERROR: unexpected aero_particle_mass dimids"
00133            stop 1
00134         end if
00135         allocate(aero_particle_mass(n_aero_particle, n_aero_species))
00136         call nc_check_msg(nf90_get_var(ncid, varid_aero_particle_mass, &
00137              aero_particle_mass), "getting variable 'aero_particle_mass'")
00138         
00139         ! read aero_comp_vol
00140         call nc_check_msg(nf90_inq_varid(ncid, "aero_comp_vol", &
00141              varid_aero_comp_vol), "getting variable ID for 'aero_comp_vol'")
00142         call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_comp_vol, &
00143              tmp_str, xtype, ndims, dimids, nAtts), &
00144              "inquiring variable 'aero_comp_vol'")
00145         if ((ndims /= 1) &
00146              .or. (dimids(1) /= dimid_aero_particle)) then
00147            write(0,*) "ERROR: unexpected aero_comp_vol dimids"
00148            stop 1
00149         end if
00150         allocate(aero_comp_vol(n_aero_particle))
00151         call nc_check_msg(nf90_get_var(ncid, varid_aero_comp_vol, &
00152              aero_comp_vol), "getting variable 'aero_comp_vol'")
00153         
00154         call nc_check_msg(nf90_close(ncid), &
00155              "closing file " // trim(in_filename))
00156         
00157         ! compute per-species total mass
00158         aero_conc = 0d0
00159         do i_part = 1,n_aero_particle
00160            do i_spec = 1,n_aero_species
00161               aero_conc(i_spec) = aero_conc(i_spec) &
00162                    + aero_particle_mass(i_part, i_spec) &
00163                    / aero_comp_vol(i_part)
00164            end do
00165         end do
00166 
00167         deallocate(aero_particle_mass)
00168         deallocate(aero_comp_vol)
00169      end if
00170 
00171      ! output data
00172      write(out_unit, '(e30.15e3)', advance='no') time
00173      do i_spec = 1,n_aero_species
00174         write(out_unit, '(e30.15e3)', advance='no') &
00175              aero_conc(i_spec)
00176      end do
00177      write(out_unit, '(a)') ''
00178 
00179   end do
00180 
00181   if (n_time == 0) then
00182      write(*,'(a,a)') 'ERROR: no input file found matching: ', &
00183           trim(in_filename)
00184      stop 1
00185   end if
00186 
00187   close(out_unit)
00188   deallocate(aero_conc)
00189 
00190 contains
00191 
00192 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00193 
00194   !> Check return status of NetCDF function calls.
00195   subroutine nc_check_msg(status, error_msg)
00196 
00197     !> Status return value.
00198     integer, intent(in) :: status
00199     !> Error message in case of failure.
00200     character(len=*), intent(in) :: error_msg
00201 
00202     if (status /= NF90_NOERR) then
00203        write(0,*) trim(error_msg) // " : " // trim(nf90_strerror(status))
00204        stop 1
00205     end if
00206 
00207   end subroutine nc_check_msg
00208 
00209 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00210 
00211 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT
00212   integer function command_argument_count()
00213     command_argument_count = iargc()
00214   end function command_argument_count
00215   subroutine get_command_argument(i, arg)
00216     integer, intent(in) :: i
00217     character(len=*), intent(out) :: arg
00218     call getarg(i, arg)
00219   end subroutine get_command_argument
00220 #endif
00221   
00222 end program extract_aero_species