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