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