PartMC 2.1.2
extract_gas.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_gas program.
00007 
00008 !> Read NetCDF output files and write out the gas mixing ratios in text
00009 !> format.
00010 program extract_gas
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_gas_species
00020   integer :: varid_time, varid_gas_species
00021   integer :: varid_gas_mixing_ratio
00022   integer :: n_gas_species
00023   character(len=1000) :: tmp_str, gas_species_names, remaining_species
00024   real(kind=dp) :: time
00025   real(kind=dp), allocatable :: gas_mixing_ratio(:)
00026   integer :: xtype, ndims, nAtts
00027   integer, dimension(nf90_max_var_dims) :: dimids
00028   integer :: ios, i_time, i_spec, status, n_time, i
00029   character(len=36) :: uuid, run_uuid
00030 
00031   ! process commandline arguments
00032   if (command_argument_count() .ne. 2) then
00033      write(6,*) 'Usage: extract_gas <netcdf_state_prefix> <output_filename>'
00034      stop 2
00035   endif
00036   call get_command_argument(1, in_prefix)
00037   call get_command_argument(2, out_filename)
00038 
00039   ! open output file
00040   open(unit=out_unit, file=out_filename, status='replace', iostat=ios)
00041   if (ios /= 0) then
00042      write(0,'(a,a,a,i4)') 'ERROR: unable to open file ', &
00043           trim(out_filename), ' for writing: ', ios
00044      stop 1
00045   end if
00046 
00047   ! write information
00048   write(*,'(a,a)') "Output file: ", trim(out_filename)
00049   write(*,'(a)') "  Each row of output is one time."
00050   write(*,'(a)') "  The columns of output are:"
00051   write(*,'(a)') "    column  1: time (s)"
00052 
00053   ! read NetCDF files
00054   i_time = 0
00055   n_time = 0
00056   do while (.true.)
00057      i_time = i_time + 1
00058      write(in_filename,'(a,i8.8,a)') trim(in_prefix), i_time, ".nc"
00059      status = nf90_open(in_filename, NF90_NOWRITE, ncid)
00060      if (status /= NF90_NOERR) then
00061         exit
00062      end if
00063      n_time = i_time
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 gas_species
00084      call nc_check_msg(nf90_inq_dimid(ncid, "gas_species", &
00085           dimid_gas_species), "getting dimension ID for 'gas_species'")
00086      call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_gas_species, &
00087           tmp_str, n_gas_species), "inquiring dimension 'gas_species'")
00088      call nc_check_msg(nf90_inq_varid(ncid, "gas_species", &
00089           varid_gas_species), "getting variable ID for 'gas_species'")
00090      call nc_check_msg(nf90_get_att(ncid, varid_gas_species, &
00091           "names", gas_species_names), &
00092           "getting attribute 'names' for variable 'gas_species'")
00093      if (i_time == 1) then
00094         remaining_species = gas_species_names
00095         do i_spec = 1,n_gas_species
00096            if (i_spec < n_gas_species) then
00097               i = index(remaining_species, ',')
00098            else
00099               i = index(remaining_species, ' ')
00100            end if
00101            write(*,'(a,i2,a,a,a)') '    column ', i_spec + 1, ': ', &
00102                 remaining_species(1:(i-1)), ' (ppb)'
00103            remaining_species = remaining_species((i+1):)
00104         end do
00105      end if
00106      
00107      ! read gas_mixing_ratio
00108      call nc_check_msg(nf90_inq_varid(ncid, "gas_mixing_ratio", &
00109           varid_gas_mixing_ratio), &
00110           "getting variable ID for 'gas_mixing_ratio'")
00111      call nc_check_msg(nf90_Inquire_Variable(ncid, varid_gas_mixing_ratio, &
00112           tmp_str, xtype, ndims, dimids, nAtts), &
00113           "inquiring variable 'gas_mixing_ratio'")
00114      if ((ndims /= 1) &
00115           .or. (dimids(1) /= dimid_gas_species)) then
00116         write(0,*) "ERROR: unexpected gas_mixing_ratio dimids"
00117         stop 1
00118      end if
00119      allocate(gas_mixing_ratio(n_gas_species))
00120      call nc_check_msg(nf90_get_var(ncid, varid_gas_mixing_ratio, &
00121           gas_mixing_ratio), "getting variable 'gas_mixing_ratio'")
00122      
00123      call nc_check_msg(nf90_close(ncid), &
00124           "closing file " // trim(in_filename))
00125 
00126      ! output data
00127      write(out_unit, '(e30.15e3)', advance='no') time
00128      do i_spec = 1,n_gas_species
00129         write(out_unit, '(e30.15e3)', advance='no') &
00130              gas_mixing_ratio(i_spec)
00131      end do
00132      write(out_unit, '(a)') ''
00133 
00134      deallocate(gas_mixing_ratio)
00135   end do
00136 
00137   if (n_time == 0) then
00138      write(*,'(a,a)') 'ERROR: no input file found matching: ', &
00139           trim(in_filename)
00140      stop 1
00141   end if
00142 
00143   close(out_unit)
00144 
00145 contains
00146 
00147 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00148 
00149   !> Check return status of NetCDF function calls.
00150   subroutine nc_check_msg(status, error_msg)
00151 
00152     !> Status return value.
00153     integer, intent(in) :: status
00154     !> Error message in case of failure.
00155     character(len=*), intent(in) :: error_msg
00156 
00157     if (status /= NF90_NOERR) then
00158        write(0,*) trim(error_msg) // " : " // trim(nf90_strerror(status))
00159        stop 1
00160     end if
00161 
00162   end subroutine nc_check_msg
00163 
00164 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00165 
00166 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT
00167   integer function command_argument_count()
00168     command_argument_count = iargc()
00169   end function command_argument_count
00170   subroutine get_command_argument(i, arg)
00171     integer, intent(in) :: i
00172     character(len=*), intent(out) :: arg
00173     call getarg(i, arg)
00174   end subroutine get_command_argument
00175 #endif
00176   
00177 end program extract_gas