PartMC
2.2.1
|
00001 ! Copyright (C) 2010-2011 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 pmc_aero_weight module. 00007 00008 !> The aero_weight_t structure and associated subroutines. 00009 module pmc_aero_weight 00010 00011 use pmc_util 00012 use pmc_constants 00013 use pmc_rand 00014 use pmc_spec_file 00015 use pmc_aero_particle 00016 use pmc_netcdf 00017 use pmc_mpi 00018 #ifdef PMC_USE_MPI 00019 use mpi 00020 #endif 00021 00022 !> Type code for an undefined or invalid weighting. 00023 integer, parameter :: AERO_WEIGHT_TYPE_INVALID = 0 00024 !> Type code for no (or flat) weighting. 00025 integer, parameter :: AERO_WEIGHT_TYPE_NONE = 1 00026 !> Type code for power function weighting. 00027 integer, parameter :: AERO_WEIGHT_TYPE_POWER = 2 00028 !> Type code for MFA weighting. 00029 integer, parameter :: AERO_WEIGHT_TYPE_MFA = 3 00030 00031 !> An aerosol size distribution weighting function. 00032 type aero_weight_t 00033 !> Computational volume (m^3). 00034 real(kind=dp) :: comp_vol 00035 !> Weight type (given by module constants). 00036 integer :: type 00037 !> Reference radius at which the weight is 1 (hard-coded at present). 00038 real(kind=dp) :: ref_radius 00039 !> Exponent for "power" weight. 00040 real(kind=dp) :: exponent 00041 end type aero_weight_t 00042 00043 contains 00044 00045 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00046 00047 !> Allocates an aero_weight. 00048 elemental subroutine aero_weight_allocate(aero_weight) 00049 00050 !> Aerosol weight. 00051 type(aero_weight_t), intent(out) :: aero_weight 00052 00053 call aero_weight_zero(aero_weight) 00054 00055 end subroutine aero_weight_allocate 00056 00057 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00058 00059 !> Allocates an aero_weight of the given size. 00060 elemental subroutine aero_weight_allocate_size(aero_weight) 00061 00062 !> Aerosol weight. 00063 type(aero_weight_t), intent(out) :: aero_weight 00064 00065 call aero_weight_zero(aero_weight) 00066 00067 end subroutine aero_weight_allocate_size 00068 00069 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00070 00071 !> Free all storage. 00072 elemental subroutine aero_weight_deallocate(aero_weight) 00073 00074 !> Aerosol weight. 00075 type(aero_weight_t), intent(in) :: aero_weight 00076 00077 end subroutine aero_weight_deallocate 00078 00079 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00080 00081 !> Zeros the contents of the \c aero_weight. 00082 elemental subroutine aero_weight_zero(aero_weight) 00083 00084 !> Aerosol weight. 00085 type(aero_weight_t), intent(inout) :: aero_weight 00086 00087 aero_weight%comp_vol = 0d0 00088 aero_weight%type = AERO_WEIGHT_TYPE_INVALID 00089 aero_weight%ref_radius = 0d0 00090 aero_weight%exponent = 0d0 00091 00092 end subroutine aero_weight_zero 00093 00094 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00095 00096 !> Zeros the \c aero_weight computational volume. 00097 elemental subroutine aero_weight_zero_comp_vol(aero_weight) 00098 00099 !> Aerosol weight. 00100 type(aero_weight_t), intent(inout) :: aero_weight 00101 00102 aero_weight%comp_vol = 0d0 00103 00104 end subroutine aero_weight_zero_comp_vol 00105 00106 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00107 00108 !> Copy an aero_weight. 00109 elemental subroutine aero_weight_copy(aero_weight_from, aero_weight_to) 00110 00111 !> Aerosol weight original. 00112 type(aero_weight_t), intent(in) :: aero_weight_from 00113 !> Aerosol weight copy. 00114 type(aero_weight_t), intent(inout) :: aero_weight_to 00115 00116 aero_weight_to%comp_vol = aero_weight_from%comp_vol 00117 aero_weight_to%type = aero_weight_from%type 00118 aero_weight_to%ref_radius = aero_weight_from%ref_radius 00119 aero_weight_to%exponent = aero_weight_from%exponent 00120 00121 end subroutine aero_weight_copy 00122 00123 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00124 00125 !> Copy an aero_weight array. 00126 subroutine aero_weight_array_copy(aero_weight_array_from, & 00127 aero_weight_array_to) 00128 00129 !> Aerosol weight array original. 00130 type(aero_weight_t), intent(in) :: aero_weight_array_from(:) 00131 !> Aerosol weight array copy. 00132 type(aero_weight_t), allocatable, intent(inout) :: aero_weight_array_to(:) 00133 00134 if (size(aero_weight_array_to) /= size(aero_weight_array_from)) then 00135 call aero_weight_deallocate(aero_weight_array_to) 00136 deallocate(aero_weight_array_to) 00137 allocate(aero_weight_array_to(size(aero_weight_array_from))) 00138 call aero_weight_allocate(aero_weight_array_to) 00139 end if 00140 call aero_weight_copy(aero_weight_array_from, aero_weight_array_to) 00141 00142 end subroutine aero_weight_array_copy 00143 00144 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00145 00146 !> Scale the computational volume by the given fraction, so 00147 !> <tt>new_comp_vol = old_comp_vol * fraction</tt>. 00148 elemental subroutine aero_weight_scale_comp_vol(aero_weight, fraction) 00149 00150 !> Aerosol weight to halve. 00151 type(aero_weight_t), intent(inout) :: aero_weight 00152 !> Fraction to scale computational volume by. 00153 real(kind=dp), intent(in) :: fraction 00154 00155 aero_weight%comp_vol = aero_weight%comp_vol * fraction 00156 00157 end subroutine aero_weight_scale_comp_vol 00158 00159 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00160 00161 !> Add the computational volume of \c aero_weight_delta to \c aero_weight. 00162 elemental subroutine aero_weight_add_comp_vol(aero_weight, aero_weight_delta) 00163 00164 !> Aerosol weight to add volume to. 00165 type(aero_weight_t), intent(inout) :: aero_weight 00166 !> Aerosol weight to add volume from. 00167 type(aero_weight_t), intent(in) :: aero_weight_delta 00168 00169 aero_weight%comp_vol = aero_weight%comp_vol + aero_weight_delta%comp_vol 00170 00171 end subroutine aero_weight_add_comp_vol 00172 00173 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00174 00175 !> Transfer the computational volume from \c aero_weight_from to \c 00176 !> aero_weight_to, weighted by \c sample_prop. 00177 elemental subroutine aero_weight_transfer_comp_vol(aero_weight_from, & 00178 aero_weight_to, sample_prop) 00179 00180 !> Aerosol weight to take volume from. 00181 type(aero_weight_t), intent(inout) :: aero_weight_from 00182 !> Aerosol weight to add volume to. 00183 type(aero_weight_t), intent(inout) :: aero_weight_to 00184 !> Proportion of from volume to transfer. 00185 real(kind=dp), intent(in) :: sample_prop 00186 00187 real(kind=dp) :: transfer_comp_vol 00188 00189 transfer_comp_vol = sample_prop * aero_weight_from%comp_vol 00190 aero_weight_to%comp_vol = aero_weight_to%comp_vol + transfer_comp_vol 00191 aero_weight_from%comp_vol = aero_weight_from%comp_vol - transfer_comp_vol 00192 00193 end subroutine aero_weight_transfer_comp_vol 00194 00195 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00196 00197 !> Compute the number concentration at a given radius (m^{-3}). 00198 real(kind=dp) function aero_weight_num_conc_at_radius(aero_weight, radius) 00199 00200 !> Aerosol weight. 00201 type(aero_weight_t), intent(in) :: aero_weight 00202 !> Radius to compute number concentration at (m). 00203 real(kind=dp), intent(in) :: radius 00204 00205 if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then 00206 aero_weight_num_conc_at_radius = 1d0 00207 elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) & 00208 .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then 00209 aero_weight_num_conc_at_radius & 00210 = (radius / aero_weight%ref_radius)**aero_weight%exponent 00211 else 00212 call die_msg(700421478, "unknown aero_weight type: " & 00213 // trim(integer_to_string(aero_weight%type))) 00214 end if 00215 aero_weight_num_conc_at_radius & 00216 = aero_weight_num_conc_at_radius / aero_weight%comp_vol 00217 00218 end function aero_weight_num_conc_at_radius 00219 00220 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00221 00222 !> Compute the number concentration for a particle (m^{-3}). 00223 real(kind=dp) function aero_weight_num_conc(aero_weight, & 00224 aero_particle) 00225 00226 !> Aerosol weight. 00227 type(aero_weight_t), intent(in) :: aero_weight 00228 !> Aerosol particle to compute number concentration for. 00229 type(aero_particle_t), intent(in) :: aero_particle 00230 00231 aero_weight_num_conc = aero_weight_num_conc_at_radius(aero_weight, & 00232 aero_particle_radius(aero_particle)) 00233 00234 end function aero_weight_num_conc 00235 00236 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00237 00238 !> Compute the number concentration for a particle (m^{-3}). 00239 real(kind=dp) function aero_weight_array_single_num_conc(aero_weight_array, & 00240 aero_particle) 00241 00242 !> Aerosol weight array. 00243 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00244 !> Aerosol particle to compute number concentration for. 00245 type(aero_particle_t), intent(in) :: aero_particle 00246 00247 aero_weight_array_single_num_conc = aero_weight_num_conc( & 00248 aero_weight_array(aero_particle%weight_group), aero_particle) 00249 00250 end function aero_weight_array_single_num_conc 00251 00252 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00253 00254 !> Compute the total number concentration at a given radius (m^3). 00255 real(kind=dp) function aero_weight_array_num_conc_at_radius( & 00256 aero_weight_array, radius) 00257 00258 !> Aerosol weight array. 00259 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00260 !> Radius to compute number concentration at (m). 00261 real(kind=dp), intent(in) :: radius 00262 00263 integer :: i_group 00264 real(kind=dp) :: num_conc(size(aero_weight_array)) 00265 00266 do i_group = 1,size(aero_weight_array) 00267 num_conc(i_group) & 00268 = aero_weight_num_conc_at_radius(aero_weight_array(i_group), & 00269 radius) 00270 end do 00271 ! harmonic mean (same as summing the computational volumes) 00272 aero_weight_array_num_conc_at_radius = 1d0 / sum(1d0 / num_conc) 00273 00274 end function aero_weight_array_num_conc_at_radius 00275 00276 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00277 00278 !> Compute the number concentration for a particle (m^{-3}). 00279 real(kind=dp) function aero_weight_array_num_conc(aero_weight_array, & 00280 aero_particle) 00281 00282 !> Aerosol weight array. 00283 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00284 !> Aerosol particle to compute number concentration for. 00285 type(aero_particle_t), intent(in) :: aero_particle 00286 00287 aero_weight_array_num_conc = aero_weight_array_num_conc_at_radius( & 00288 aero_weight_array, aero_particle_radius(aero_particle)) 00289 00290 end function aero_weight_array_num_conc 00291 00292 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00293 00294 !> Compute the radius at a given number concentration (m). Inverse 00295 !> of aero_weight_num_conc_at_radius(). 00296 real(kind=dp) function aero_weight_radius_at_num_conc(aero_weight, num_conc) 00297 00298 !> Aerosol weight. 00299 type(aero_weight_t), intent(in) :: aero_weight 00300 !> Number concentration to compute the radius at (m^{-3}). 00301 real(kind=dp), intent(in) :: num_conc 00302 00303 if (aero_weight%type == AERO_WEIGHT_TYPE_NONE) then 00304 call die_msg(545688537, "cannot invert weight type 'none'") 00305 elseif ((aero_weight%type == AERO_WEIGHT_TYPE_POWER) & 00306 .or. (aero_weight%type == AERO_WEIGHT_TYPE_MFA)) then 00307 call assert_msg(902242996, aero_weight%exponent /= 0d0, & 00308 "cannot invert weight with zero exponent") 00309 aero_weight_radius_at_num_conc = exp(log(aero_weight%ref_radius) & 00310 + log(num_conc * aero_weight%comp_vol) / aero_weight%exponent) 00311 else 00312 call die_msg(307766845, "unknown aero_weight type: " & 00313 // trim(integer_to_string(aero_weight%type))) 00314 end if 00315 00316 end function aero_weight_radius_at_num_conc 00317 00318 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00319 00320 !> Check whether a given aero_weight array is flat in total. 00321 logical function aero_weight_array_check_flat(aero_weight_array) 00322 00323 !> Aerosol weight array. 00324 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00325 00326 integer :: i_group 00327 00328 ! check we know about all the weight types 00329 do i_group = 1,size(aero_weight_array) 00330 call assert(269952052, & 00331 (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_NONE) & 00332 .or. (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_POWER) & 00333 .or. (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_MFA)) 00334 if (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_NONE) then 00335 call assert(853998284, aero_weight_array(i_group)%exponent == 0d0) 00336 end if 00337 if (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_MFA) then 00338 call assert(829651126, aero_weight_array(i_group)%exponent == -3d0) 00339 end if 00340 end do 00341 00342 if (abs(sum(aero_weight_array%exponent)) & 00343 < 1d-20 * sum(abs(aero_weight_array%exponent))) then 00344 aero_weight_array_check_flat = .true. 00345 else 00346 aero_weight_array_check_flat = .false. 00347 end if 00348 00349 end function aero_weight_array_check_flat 00350 00351 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00352 00353 !> Determine whether all weight functions in an array are monotone 00354 !> increasing, monotone decreasing, or neither. 00355 subroutine aero_weight_array_check_monotonicity(aero_weight_array, & 00356 monotone_increasing, monotone_decreasing) 00357 00358 !> Aerosol weight array. 00359 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00360 !> Whether all weights are monotone increasing. 00361 logical, intent(out) :: monotone_increasing 00362 !> Whether all weights are monotone decreasing. 00363 logical, intent(out) :: monotone_decreasing 00364 00365 integer :: i_group 00366 00367 monotone_increasing = .true. 00368 monotone_decreasing = .true. 00369 do i_group = 1,size(aero_weight_array) 00370 ! check we know about all the weight types 00371 call assert(610698264, & 00372 (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_NONE) & 00373 .or. (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_POWER) & 00374 .or. (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_MFA)) 00375 if (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_POWER) then 00376 if (aero_weight_array(i_group)%exponent < 0d0) then 00377 monotone_increasing = .false. 00378 end if 00379 if (aero_weight_array(i_group)%exponent > 0d0) then 00380 monotone_decreasing = .false. 00381 end if 00382 end if 00383 if (aero_weight_array(i_group)%type == AERO_WEIGHT_TYPE_MFA) then 00384 monotone_increasing = .false. 00385 end if 00386 end do 00387 00388 end subroutine aero_weight_array_check_monotonicity 00389 00390 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00391 00392 !> Compute the maximum and minimum number concentrations between the 00393 !> given radii. 00394 subroutine aero_weight_array_minmax_num_conc(aero_weight_array, radius_1, & 00395 radius_2, num_conc_min, num_conc_max) 00396 00397 !> Aerosol weight. 00398 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00399 !> First radius. 00400 real(kind=dp), intent(in) :: radius_1 00401 !> Second radius. 00402 real(kind=dp), intent(in) :: radius_2 00403 !> Minimum number concentration. 00404 real(kind=dp), intent(out) :: num_conc_min 00405 !> Maximum number concentration. 00406 real(kind=dp), intent(out) :: num_conc_max 00407 00408 real(kind=dp) :: num_conc_1, num_conc_2 00409 logical :: monotone_increasing, monotone_decreasing 00410 00411 call aero_weight_array_check_monotonicity(aero_weight_array, & 00412 monotone_increasing, monotone_decreasing) 00413 call assert(857727714, monotone_increasing .or. monotone_decreasing) 00414 00415 num_conc_1 = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00416 radius_1) 00417 num_conc_2 = aero_weight_array_num_conc_at_radius(aero_weight_array, & 00418 radius_2) 00419 num_conc_min = min(num_conc_1, num_conc_2) 00420 num_conc_max = max(num_conc_1, num_conc_2) 00421 00422 end subroutine aero_weight_array_minmax_num_conc 00423 00424 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00425 00426 !> Choose a random group at the given radius, with probability 00427 !> proportional to group volume at that radius. 00428 integer function aero_weight_array_rand_group(aero_weight_array, radius) 00429 00430 !> Aerosol weight array. 00431 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00432 !> Radius to sample group at (m). 00433 real(kind=dp), intent(in) :: radius 00434 00435 real(kind=dp) :: comp_vols(size(aero_weight_array)) 00436 integer :: i 00437 00438 do i = 1,size(aero_weight_array) 00439 comp_vols(i) & 00440 = 1d0 / aero_weight_num_conc_at_radius(aero_weight_array(i), & 00441 radius) 00442 end do 00443 aero_weight_array_rand_group = sample_cts_pdf(comp_vols) 00444 00445 end function aero_weight_array_rand_group 00446 00447 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00448 00449 !> Read an aero_weight from a spec file. 00450 subroutine spec_file_read_aero_weight(file, aero_weight) 00451 00452 !> Spec file. 00453 type(spec_file_t), intent(inout) :: file 00454 !> Aerosol weight. 00455 type(aero_weight_t), intent(inout) :: aero_weight 00456 00457 character(len=SPEC_LINE_MAX_VAR_LEN) :: weight_type 00458 00459 !> \page input_format_aero_weight Input File Format: Aerosol Weighting Function 00460 !! 00461 !! For efficiency the aerosol population can be weighted so that 00462 !! the true number distribution \f$n(D)\f$ is given by 00463 !! \f[ n(D) = w(D) c(D) \f] 00464 !! where \f$w(D)\f$ is a fixed weighting function, \f$c(D)\f$ is 00465 !! the computational (simulated) number distribution, and \f$D\f$ 00466 !! is the diameter. Thus a large value of \f$w(D)\f$ means that 00467 !! relatively few computational particles are used at diameter 00468 !! \f$D\f$, while a small value of \f$w(D)\f$ means that 00469 !! relatively many computational particles will be used at that 00470 !! diameter. 00471 !! 00472 !! The aerosol weighting function is specified by the parameters: 00473 !! - \b weight (string): the type of weighting function --- must 00474 !! be one of: "none" for no weighting (\f$w(D) = 1\f$); 00475 !! "power" for a power-law weighting (\f$w(D) = D^\alpha\f$), 00476 !! or "mfa" for the mass flow algorithm weighting (\f$w(D) = 00477 !! D^{-3}\f$ with dependent coagulation particle removal) 00478 !! - if the \c weight is \c power then the next parameter is: 00479 !! - \b exponent (real, dimensionless): the exponent 00480 !! \f$\alpha\f$ in the power law relationship --- setting 00481 !! the \c exponent to 0 is equivalent to no weighting, while 00482 !! setting the exponent negative uses more computational 00483 !! particles at larger diameters and setting the exponent 00484 !! positive uses more computational particles at smaller 00485 !! diameters; in practice exponents between 0 and -3 are 00486 !! most useful 00487 !! 00488 !! See also: 00489 !! - \ref spec_file_format --- the input file text format 00490 00491 call spec_file_read_string(file, 'weight', weight_type) 00492 if (trim(weight_type) == 'none') then 00493 aero_weight%type = AERO_WEIGHT_TYPE_NONE 00494 elseif (trim(weight_type) == 'power') then 00495 aero_weight%type = AERO_WEIGHT_TYPE_POWER 00496 aero_weight%ref_radius = 1d0 00497 call spec_file_read_real(file, 'exponent', aero_weight%exponent) 00498 elseif (trim(weight_type) == 'mfa') then 00499 aero_weight%type = AERO_WEIGHT_TYPE_MFA 00500 aero_weight%ref_radius = 1d0 00501 aero_weight%exponent = -3d0 00502 else 00503 call spec_file_die_msg(456342050, file, "unknown weight_type: " & 00504 // trim(weight_type)) 00505 end if 00506 00507 end subroutine spec_file_read_aero_weight 00508 00509 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00510 00511 !> Determines the number of bytes required to pack the given value. 00512 integer function pmc_mpi_pack_size_aero_weight(val) 00513 00514 !> Value to pack. 00515 type(aero_weight_t), intent(in) :: val 00516 00517 pmc_mpi_pack_size_aero_weight = & 00518 pmc_mpi_pack_size_real(val%comp_vol) & 00519 + pmc_mpi_pack_size_integer(val%type) & 00520 + pmc_mpi_pack_size_real(val%ref_radius) & 00521 + pmc_mpi_pack_size_real(val%exponent) 00522 00523 end function pmc_mpi_pack_size_aero_weight 00524 00525 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00526 00527 !> Packs the given value into the buffer, advancing position. 00528 subroutine pmc_mpi_pack_aero_weight(buffer, position, val) 00529 00530 !> Memory buffer. 00531 character, intent(inout) :: buffer(:) 00532 !> Current buffer position. 00533 integer, intent(inout) :: position 00534 !> Value to pack. 00535 type(aero_weight_t), intent(in) :: val 00536 00537 #ifdef PMC_USE_MPI 00538 integer :: prev_position 00539 00540 prev_position = position 00541 call pmc_mpi_pack_real(buffer, position, val%comp_vol) 00542 call pmc_mpi_pack_integer(buffer, position, val%type) 00543 call pmc_mpi_pack_real(buffer, position, val%ref_radius) 00544 call pmc_mpi_pack_real(buffer, position, val%exponent) 00545 call assert(579699255, & 00546 position - prev_position <= pmc_mpi_pack_size_aero_weight(val)) 00547 #endif 00548 00549 end subroutine pmc_mpi_pack_aero_weight 00550 00551 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00552 00553 !> Unpacks the given value from the buffer, advancing position. 00554 subroutine pmc_mpi_unpack_aero_weight(buffer, position, val) 00555 00556 !> Memory buffer. 00557 character, intent(inout) :: buffer(:) 00558 !> Current buffer position. 00559 integer, intent(inout) :: position 00560 !> Value to pack. 00561 type(aero_weight_t), intent(inout) :: val 00562 00563 #ifdef PMC_USE_MPI 00564 integer :: prev_position 00565 00566 prev_position = position 00567 call pmc_mpi_unpack_real(buffer, position, val%comp_vol) 00568 call pmc_mpi_unpack_integer(buffer, position, val%type) 00569 call pmc_mpi_unpack_real(buffer, position, val%ref_radius) 00570 call pmc_mpi_unpack_real(buffer, position, val%exponent) 00571 call assert(639056899, & 00572 position - prev_position <= pmc_mpi_pack_size_aero_weight(val)) 00573 #endif 00574 00575 end subroutine pmc_mpi_unpack_aero_weight 00576 00577 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00578 00579 !> Write the \c aero_weight dimension to the given NetCDF file if it 00580 !> is not already present and in any case return the associated 00581 !> dimid. 00582 subroutine aero_weight_netcdf_dim_aero_weight(aero_weight_array, ncid, & 00583 dimid_aero_weight) 00584 00585 !> Aero_weight structure array. 00586 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00587 !> NetCDF file ID, in data mode. 00588 integer, intent(in) :: ncid 00589 !> Dimid of the species dimension. 00590 integer, intent(out) :: dimid_aero_weight 00591 00592 integer :: status, i_weight, n_weight 00593 integer :: varid_aero_weight 00594 integer :: aero_weight_centers(size(aero_weight_array)) 00595 00596 ! try to get the dimension ID 00597 status = nf90_inq_dimid(ncid, "aero_weight", dimid_aero_weight) 00598 if (status == NF90_NOERR) return 00599 if (status /= NF90_EBADDIM) call pmc_nc_check(status) 00600 00601 ! dimension not defined, so define now define it 00602 call pmc_nc_check(nf90_redef(ncid)) 00603 00604 n_weight = size(aero_weight_array) 00605 00606 call pmc_nc_check(nf90_def_dim(ncid, "aero_weight", n_weight, & 00607 dimid_aero_weight)) 00608 call pmc_nc_check(nf90_def_var(ncid, "aero_weight", NF90_INT, & 00609 dimid_aero_weight, varid_aero_weight)) 00610 call pmc_nc_check(nf90_put_att(ncid, varid_aero_weight, "description", & 00611 "dummy dimension variable (no useful value)")) 00612 00613 call pmc_nc_check(nf90_enddef(ncid)) 00614 00615 do i_weight = 1,n_weight 00616 aero_weight_centers(i_weight) = i_weight 00617 end do 00618 call pmc_nc_check(nf90_put_var(ncid, varid_aero_weight, & 00619 aero_weight_centers)) 00620 00621 end subroutine aero_weight_netcdf_dim_aero_weight 00622 00623 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00624 00625 !> Write full state. 00626 subroutine aero_weight_array_output_netcdf(aero_weight_array, ncid) 00627 00628 !> Aero weight to write. 00629 type(aero_weight_t), intent(in) :: aero_weight_array(:) 00630 !> NetCDF file ID, in data mode. 00631 integer, intent(in) :: ncid 00632 00633 integer :: dimid_aero_weight 00634 00635 !> \page output_format_aero_weight_array Output File Format: Aerosol Weighting Functions 00636 !! 00637 !! The aerosol weighting function NetCDF dimensions are: 00638 !! - \b aero_weight: number of aerosol weighting functions 00639 !! 00640 !! The aerosol weighting function NetCDF variables are: 00641 !! - \b weight_comp_vol (unit m^3, dim \c aero_weight): the 00642 !! computational volume associated with each weighting 00643 !! function 00644 !! - \b weight_type (no unit, dim \c aero_weight): the type of 00645 !! each weighting function, with 0 = invalid weight, 1 = no 00646 !! weight (\f$w(D) = 1\f$), 2 = power weight (\f$w(D) = 00647 !! (D/D_0)^\alpha\f$), 3 = MFA weight (\f$w(D) = 00648 !! (D/D_0)^{-3}\f$) 00649 !! - \b weight_exponent (no unit, dim \c aero_weight): for each 00650 !! weighting function, specifies the exponent \f$\alpha\f$ for 00651 !! the power \c weight_type, the value -3 for the MFA \c 00652 !! weight_type, and zero for any other \c weight_type 00653 00654 call aero_weight_netcdf_dim_aero_weight(aero_weight_array, ncid, & 00655 dimid_aero_weight) 00656 00657 call pmc_nc_write_real_1d(ncid, aero_weight_array%comp_vol, & 00658 "weight_comp_vol", (/ dimid_aero_weight /), unit="m^3", & 00659 description="computational volume for each weighting function") 00660 call pmc_nc_write_integer_1d(ncid, aero_weight_array%type, "weight_type", & 00661 (/ dimid_aero_weight /), & 00662 description="type of each aerosol weighting function: 0 = invalid, " & 00663 // "1 = none (w(D) = 1), 2 = power (w(D) = (D/D_0)^alpha), " & 00664 // "3 = MFA (mass flow) (w(D) = (D/D_0)^(-3))") 00665 call pmc_nc_write_real_1d(ncid, aero_weight_array%exponent, & 00666 "weight_exponent", (/ dimid_aero_weight /), unit="1", & 00667 description="exponent alpha for the power weight_type, " & 00668 // "set to -3 for MFA, and zero otherwise") 00669 00670 end subroutine aero_weight_array_output_netcdf 00671 00672 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00673 00674 !> Read full state. 00675 subroutine aero_weight_array_input_netcdf(aero_weight_array, ncid) 00676 00677 !> Environment state to read. 00678 type(aero_weight_t), intent(inout), allocatable :: aero_weight_array(:) 00679 !> NetCDF file ID, in data mode. 00680 integer, intent(in) :: ncid 00681 00682 integer :: dimid_aero_weight, n_weight 00683 character(len=1000) :: name 00684 real(kind=dp), allocatable :: comp_vol(:), exponent(:) 00685 integer, allocatable :: type(:) 00686 00687 call pmc_nc_check(nf90_inq_dimid(ncid, "aero_weight", & 00688 dimid_aero_weight)) 00689 call pmc_nc_check(nf90_Inquire_Dimension(ncid, & 00690 dimid_aero_weight, name, n_weight)) 00691 call assert(719221386, n_weight < 1000) 00692 00693 allocate(comp_vol(n_weight)) 00694 allocate(type(n_weight)) 00695 allocate(exponent(n_weight)) 00696 00697 call pmc_nc_read_real_1d(ncid, comp_vol, "weight_comp_vol") 00698 call pmc_nc_read_integer_1d(ncid, type, "weight_type") 00699 call pmc_nc_read_real_1d(ncid, exponent, "weight_exponent") 00700 00701 call assert(309191498, size(comp_vol) == size(type)) 00702 call assert(588649520, size(comp_vol) == size(exponent)) 00703 00704 call aero_weight_deallocate(aero_weight_array) 00705 deallocate(aero_weight_array) 00706 allocate(aero_weight_array(size(comp_vol))) 00707 call aero_weight_allocate(aero_weight_array) 00708 00709 aero_weight_array%comp_vol = comp_vol 00710 aero_weight_array%type = type 00711 aero_weight_array%ref_radius = 1d0 00712 aero_weight_array%exponent = exponent 00713 00714 deallocate(comp_vol) 00715 deallocate(type) 00716 deallocate(exponent) 00717 00718 end subroutine aero_weight_array_input_netcdf 00719 00720 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00721 00722 end module pmc_aero_weight