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_sectional_aero_total program. 00007 00008 !> Read NetCDF sectional output files and write out the total aerosol 00009 !> number and mass concentrations in text format. 00010 program extract_sectional_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_diam 00020 integer :: varid_time, varid_aero_species 00021 integer :: varid_aero_diam, varid_aero_diam_widths 00022 integer :: varid_aero_mass_concentration, varid_aero_number_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_diam(:) 00027 real(kind=dp), allocatable :: aero_diam_widths(:) 00028 real(kind=dp), allocatable :: aero_mass_concentration(:,:) 00029 real(kind=dp), allocatable :: aero_number_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 real(kind=dp) :: num_conc, mass_conc 00035 character(len=36) :: uuid, run_uuid 00036 00037 ! process commandline arguments 00038 if (command_argument_count() .ne. 2) then 00039 write(6,*) 'Usage: extract_sectional_aero_total' & 00040 // ' <netcdf_state_prefix> <output_filename>' 00041 stop 2 00042 endif 00043 call get_command_argument(1, in_prefix) 00044 call get_command_argument(2, out_filename) 00045 00046 ! open output file 00047 open(unit=out_unit, file=out_filename, status='replace', iostat=ios) 00048 if (ios /= 0) then 00049 write(0,'(a,a,a,i4)') 'ERROR: unable to open file ', & 00050 trim(out_filename), ' for writing: ', ios 00051 stop 1 00052 end if 00053 00054 ! write information 00055 write(*,*) "Output file array A has:" 00056 write(*,*) " A(i, 1) = time(i) (s)" 00057 write(*,*) " A(i, 2) = total number concentration at time(i) (#/m^3)" 00058 write(*,*) " A(i, 3) = total mass concentration at time(i) (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 00096 ! read aero_diam dimension 00097 call nc_check_msg(nf90_inq_dimid(ncid, "aero_diam", & 00098 dimid_aero_diam), "getting dimension ID for 'aero_diam'") 00099 call nc_check_msg(nf90_Inquire_Dimension(ncid, dimid_aero_diam, & 00100 tmp_str, n_bin), "inquiring dimension 'aero_diam'") 00101 00102 ! read aero_diam variable 00103 call nc_check_msg(nf90_inq_varid(ncid, "aero_diam", & 00104 varid_aero_diam), & 00105 "getting variable ID for 'aero_diam'") 00106 call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_diam, & 00107 tmp_str, xtype, ndims, dimids, nAtts), & 00108 "inquiring variable 'aero_diam'") 00109 if ((ndims /= 1) & 00110 .or. (dimids(1) /= dimid_aero_diam)) then 00111 write(0,*) "ERROR: unexpected aero_diam dimids" 00112 stop 1 00113 end if 00114 allocate(aero_diam(n_bin)) 00115 call nc_check_msg(nf90_get_var(ncid, varid_aero_diam, & 00116 aero_diam), "getting variable 'aero_diam'") 00117 00118 ! read aero_diam_widths variable 00119 call nc_check_msg(nf90_inq_varid(ncid, "aero_diam_widths", & 00120 varid_aero_diam_widths), & 00121 "getting variable ID 'aero_diam_widths'") 00122 call nc_check_msg(nf90_Inquire_Variable(ncid, varid_aero_diam_widths, & 00123 tmp_str, xtype, ndims, dimids, nAtts), & 00124 "inquiring variable 'aero_diam_widths'") 00125 if ((ndims /= 1) & 00126 .or. (dimids(1) /= dimid_aero_diam)) then 00127 write(0,*) "ERROR: unexpected aero_diam_widths dimids" 00128 stop 1 00129 end if 00130 allocate(aero_diam_widths(n_bin)) 00131 call nc_check_msg(nf90_get_var(ncid, varid_aero_diam_widths, & 00132 aero_diam_widths), "getting variable 'aero_diam_widths'") 00133 00134 ! read aero_number_concentration 00135 call nc_check_msg(nf90_inq_varid(ncid, "aero_number_concentration", & 00136 varid_aero_number_concentration), & 00137 "getting variable ID for 'aero_number_concentration'") 00138 call nc_check_msg(nf90_Inquire_Variable(ncid, & 00139 varid_aero_number_concentration, tmp_str, xtype, ndims, & 00140 dimids, nAtts), "inquiring variable 'aero_number_concentration'") 00141 if ((ndims /= 1) & 00142 .or. (dimids(1) /= dimid_aero_diam)) then 00143 write(0,*) "ERROR: unexpected aero_number_concentration dimids" 00144 stop 1 00145 end if 00146 allocate(aero_number_concentration(n_bin)) 00147 call nc_check_msg(nf90_get_var(ncid, varid_aero_number_concentration, & 00148 aero_number_concentration), & 00149 "getting variable 'aero_number_concentration'") 00150 00151 ! read aero_mass_concentration 00152 call nc_check_msg(nf90_inq_varid(ncid, "aero_mass_concentration", & 00153 varid_aero_mass_concentration), & 00154 "getting variable ID for 'aero_mass_concentration'") 00155 call nc_check_msg(nf90_Inquire_Variable(ncid, & 00156 varid_aero_mass_concentration, tmp_str, xtype, ndims, & 00157 dimids, nAtts), "inquiring variable 'aero_mass_concentration'") 00158 if ((ndims /= 2) & 00159 .or. (dimids(1) /= dimid_aero_diam) & 00160 .or. (dimids(2) /= dimid_aero_species)) then 00161 write(0,*) "ERROR: unexpected aero_mass_concentration dimids" 00162 stop 1 00163 end if 00164 allocate(aero_mass_concentration(n_bin, n_aero_species)) 00165 call nc_check_msg(nf90_get_var(ncid, varid_aero_mass_concentration, & 00166 aero_mass_concentration), & 00167 "getting variable 'aero_mass_concentration'") 00168 00169 call nc_check_msg(nf90_close(ncid), & 00170 "closing file " // trim(in_filename)) 00171 00172 ! compute number and mass concentrations 00173 num_conc = sum(aero_number_concentration * aero_diam_widths) 00174 mass_conc = 0d0 00175 do i_spec = 1,n_aero_species 00176 mass_conc = mass_conc + sum(aero_mass_concentration(:,i_spec) & 00177 * aero_diam_widths) 00178 end do 00179 00180 deallocate(aero_diam) 00181 deallocate(aero_diam_widths) 00182 deallocate(aero_mass_concentration) 00183 deallocate(aero_number_concentration) 00184 00185 ! output data 00186 write(out_unit, '(3e30.15e3)') time, num_conc, mass_conc 00187 00188 end do 00189 00190 if (n_time == 0) then 00191 write(*,'(a,a)') 'ERROR: no input file found matching: ', & 00192 trim(in_filename) 00193 stop 1 00194 end if 00195 00196 close(out_unit) 00197 00198 contains 00199 00200 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00201 00202 !> Check return status of NetCDF function calls. 00203 subroutine nc_check_msg(status, error_msg) 00204 00205 !> Status return value. 00206 integer, intent(in) :: status 00207 !> Error message in case of failure. 00208 character(len=*), intent(in) :: error_msg 00209 00210 if (status /= NF90_NOERR) then 00211 write(0,*) trim(error_msg) // " : " // trim(nf90_strerror(status)) 00212 stop 1 00213 end if 00214 00215 end subroutine nc_check_msg 00216 00217 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00218 00219 #ifdef DEFINE_LOCAL_COMMAND_ARGUMENT 00220 integer function command_argument_count() 00221 command_argument_count = iargc() 00222 end function command_argument_count 00223 subroutine get_command_argument(i, arg) 00224 integer, intent(in) :: i 00225 character(len=*), intent(out) :: arg 00226 call getarg(i, arg) 00227 end subroutine get_command_argument 00228 #endif 00229 00230 end program extract_sectional_aero_total