PartMC 2.1.4
|
00001 ! Copyright (C) 2005-2011 Nicole Riemer and 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_env_state module. 00007 00008 !> The env_state_t structure and associated subroutines. 00009 module pmc_env_state 00010 00011 use pmc_gas_state 00012 use pmc_aero_dist 00013 use pmc_constants 00014 use pmc_aero_data 00015 use pmc_aero_weight 00016 use pmc_aero_particle 00017 use pmc_aero_binned 00018 use pmc_util 00019 use pmc_gas_data 00020 use pmc_bin_grid 00021 use pmc_aero_state 00022 use pmc_spec_file 00023 use pmc_mpi 00024 use pmc_netcdf 00025 #ifdef PMC_USE_MPI 00026 use mpi 00027 #endif 00028 00029 !> Current environment state. 00030 !! 00031 !! All quantities are instantaneous, describing the state at a 00032 !! particular instant of time. Constant data and other data not 00033 !! associated with the current environment state is store in 00034 !! env_data_t. 00035 !! 00036 !! The emissions and dilution are both described by pairs of a state 00037 !! and a rate. The product of these gives the actual emissions or 00038 !! dilution with units quantity per time. One way to think about 00039 !! this is to set the rate to 1/3600 and then regard the state as an 00040 !! amount per hour, etc. 00041 type env_state_t 00042 !> Temperature (K). 00043 real(kind=dp) :: temp 00044 !> Relative humidity (1). 00045 real(kind=dp) :: rel_humid 00046 !> Ambient pressure (Pa). 00047 real(kind=dp) :: pressure 00048 !> Longitude (degrees). 00049 real(kind=dp) :: longitude 00050 !> Latitude (degrees). 00051 real(kind=dp) :: latitude 00052 !> Altitude (m). 00053 real(kind=dp) :: altitude 00054 !> Start time (s since 00:00 UTC on \c start_day). 00055 real(kind=dp) :: start_time 00056 !> Start day of year (UTC). 00057 integer :: start_day 00058 !> Time since \c start_time (s). 00059 real(kind=dp) :: elapsed_time 00060 !> Solar zenith angle (radians from zenith). 00061 real(kind=dp) :: solar_zenith_angle 00062 !> Box height (m). 00063 real(kind=dp) :: height 00064 !> Gas emissions. 00065 type(gas_state_t) :: gas_emissions 00066 !> Gas emisssion rate (s^{-1}). 00067 real(kind=dp) :: gas_emission_rate 00068 !> Background gas mixing ratios. 00069 type(gas_state_t) :: gas_background 00070 !> Gas-background dilution rate (s^{-1}). 00071 real(kind=dp) :: gas_dilution_rate 00072 !> Aerosol emissions. 00073 type(aero_dist_t) :: aero_emissions 00074 !> Aerosol emisssion rate (s^{-1}). 00075 real(kind=dp) :: aero_emission_rate 00076 !> Aerosol background. 00077 type(aero_dist_t) :: aero_background 00078 !> Aero-background dilute rate (s^{-1}). 00079 real(kind=dp) :: aero_dilution_rate 00080 end type env_state_t 00081 00082 contains 00083 00084 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00085 00086 !> Allocate an empty environment. 00087 subroutine env_state_allocate(env_state) 00088 00089 !> Environment. 00090 type(env_state_t), intent(out) :: env_state 00091 00092 env_state%temp = 0d0 00093 env_state%rel_humid = 0d0 00094 env_state%pressure = 0d0 00095 env_state%longitude = 0d0 00096 env_state%latitude = 0d0 00097 env_state%altitude = 0d0 00098 env_state%start_time = 0d0 00099 env_state%start_day = 0 00100 env_state%elapsed_time = 0d0 00101 env_state%solar_zenith_angle = 0d0 00102 env_state%height = 0d0 00103 00104 call gas_state_allocate(env_state%gas_emissions) 00105 call gas_state_allocate(env_state%gas_background) 00106 env_state%gas_emission_rate = 0d0 00107 env_state%gas_dilution_rate = 0d0 00108 call aero_dist_allocate(env_state%aero_emissions) 00109 call aero_dist_allocate(env_state%aero_background) 00110 env_state%aero_emission_rate = 0d0 00111 env_state%aero_dilution_rate = 0d0 00112 00113 end subroutine env_state_allocate 00114 00115 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00116 00117 !> Free all storage. 00118 subroutine env_state_deallocate(env_state) 00119 00120 !> Environment. 00121 type(env_state_t), intent(inout) :: env_state 00122 00123 call gas_state_deallocate(env_state%gas_emissions) 00124 call gas_state_deallocate(env_state%gas_background) 00125 call aero_dist_deallocate(env_state%aero_emissions) 00126 call aero_dist_deallocate(env_state%aero_background) 00127 00128 end subroutine env_state_deallocate 00129 00130 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00131 00132 !> env_state += env_state_delta 00133 subroutine env_state_add(env_state, env_state_delta) 00134 00135 !> Environment. 00136 type(env_state_t), intent(inout) :: env_state 00137 !> Increment. 00138 type(env_state_t), intent(in) :: env_state_delta 00139 00140 env_state%temp = env_state%temp + env_state_delta%temp 00141 env_state%rel_humid = env_state%rel_humid + env_state_delta%rel_humid 00142 env_state%pressure = env_state%pressure + env_state_delta%pressure 00143 env_state%longitude = env_state%longitude + env_state_delta%longitude 00144 env_state%latitude = env_state%latitude + env_state_delta%latitude 00145 env_state%altitude = env_state%altitude + env_state_delta%altitude 00146 env_state%start_time = env_state%start_time + env_state_delta%start_time 00147 env_state%start_day = env_state%start_day + env_state_delta%start_day 00148 env_state%elapsed_time = env_state%elapsed_time & 00149 + env_state_delta%elapsed_time 00150 env_state%solar_zenith_angle = env_state%solar_zenith_angle & 00151 + env_state_delta%solar_zenith_angle 00152 env_state%height = env_state%height + env_state_delta%height 00153 call gas_state_add(env_state%gas_emissions, env_state_delta%gas_emissions) 00154 env_state%gas_emission_rate = env_state%gas_emission_rate & 00155 + env_state_delta%gas_emission_rate 00156 call gas_state_add(env_state%gas_background, & 00157 env_state_delta%gas_background) 00158 env_state%gas_dilution_rate = env_state%gas_dilution_rate & 00159 + env_state_delta%gas_dilution_rate 00160 env_state%aero_emission_rate = env_state%aero_emission_rate & 00161 + env_state_delta%aero_emission_rate 00162 env_state%aero_dilution_rate = env_state%aero_dilution_rate & 00163 + env_state_delta%aero_dilution_rate 00164 00165 end subroutine env_state_add 00166 00167 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00168 00169 !> env_state *= alpha 00170 subroutine env_state_scale(env_state, alpha) 00171 00172 !> Environment. 00173 type(env_state_t), intent(inout) :: env_state 00174 !> Scale factor. 00175 real(kind=dp), intent(in) :: alpha 00176 00177 env_state%temp = env_state%temp * alpha 00178 env_state%rel_humid = env_state%rel_humid * alpha 00179 env_state%pressure = env_state%pressure * alpha 00180 env_state%longitude = env_state%longitude * alpha 00181 env_state%latitude = env_state%latitude * alpha 00182 env_state%altitude = env_state%altitude * alpha 00183 env_state%start_time = env_state%start_time * alpha 00184 env_state%start_day = nint(real(env_state%start_day, kind=dp) * alpha) 00185 env_state%elapsed_time = env_state%elapsed_time * alpha 00186 env_state%solar_zenith_angle = env_state%solar_zenith_angle * alpha 00187 env_state%height = env_state%height * alpha 00188 call gas_state_scale(env_state%gas_emissions, alpha) 00189 env_state%gas_emission_rate = env_state%gas_emission_rate * alpha 00190 call gas_state_scale(env_state%gas_background, alpha) 00191 env_state%gas_dilution_rate = env_state%gas_dilution_rate * alpha 00192 00193 end subroutine env_state_scale 00194 00195 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00196 00197 !> env_to = env_from 00198 subroutine env_state_copy(env_from, env_to) 00199 00200 !> Original. 00201 type(env_state_t), intent(in) :: env_from 00202 !> Destination. 00203 type(env_state_t), intent(inout) :: env_to 00204 00205 env_to%temp = env_from%temp 00206 env_to%rel_humid = env_from%rel_humid 00207 env_to%pressure = env_from%pressure 00208 env_to%longitude = env_from%longitude 00209 env_to%latitude = env_from%latitude 00210 env_to%altitude = env_from%altitude 00211 env_to%start_time = env_from%start_time 00212 env_to%start_day = env_from%start_day 00213 env_to%elapsed_time = env_from%elapsed_time 00214 env_to%solar_zenith_angle = env_from%solar_zenith_angle 00215 env_to%height = env_from%height 00216 call gas_state_copy(env_from%gas_emissions, env_to%gas_emissions) 00217 env_to%gas_emission_rate = env_from%gas_emission_rate 00218 call gas_state_copy(env_from%gas_background, env_to%gas_background) 00219 env_to%gas_dilution_rate = env_from%gas_dilution_rate 00220 call aero_dist_copy(env_from%aero_emissions, env_to%aero_emissions) 00221 env_to%aero_emission_rate = env_from%aero_emission_rate 00222 call aero_dist_copy(env_from%aero_background, env_to%aero_background) 00223 env_to%aero_dilution_rate = env_from%aero_dilution_rate 00224 00225 end subroutine env_state_copy 00226 00227 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00228 00229 !> Adds the given water volume to the water vapor and updates all 00230 !> environment quantities. 00231 subroutine env_state_change_water_volume(env_state, aero_data, dv) 00232 00233 !> Environment state to update. 00234 type(env_state_t), intent(inout) :: env_state 00235 !> Aero_data constants. 00236 type(aero_data_t), intent(in) :: aero_data 00237 !> Volume concentration of water added (m^3/m^3). 00238 real(kind=dp), intent(in) :: dv 00239 00240 real(kind=dp) pmv ! ambient water vapor pressure (Pa) 00241 real(kind=dp) mv ! ambient water vapor density (kg m^{-3}) 00242 ! pmv and mv are related by the factor molec_weight/(R*T) 00243 real(kind=dp) dmv ! change of water density (kg m^{-3}) 00244 00245 dmv = dv * aero_data%density(aero_data%i_water) 00246 pmv = env_state_sat_vapor_pressure(env_state) * env_state%rel_humid 00247 mv = aero_data%molec_weight(aero_data%i_water) & 00248 / (const%univ_gas_const*env_state%temp) * pmv 00249 mv = mv - dmv 00250 if (mv < 0d0) then 00251 call warn_msg(980320483, "relative humidity tried to go negative") 00252 mv = 0d0 00253 end if 00254 env_state%rel_humid = const%univ_gas_const * env_state%temp & 00255 / aero_data%molec_weight(aero_data%i_water) * mv & 00256 / env_state_sat_vapor_pressure(env_state) 00257 00258 end subroutine env_state_change_water_volume 00259 00260 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00261 00262 !> Computes the current saturation vapor pressure (Pa). 00263 real(kind=dp) function env_state_sat_vapor_pressure(env_state) 00264 00265 !> Environment state. 00266 type(env_state_t), intent(in) :: env_state 00267 00268 env_state_sat_vapor_pressure = const%water_eq_vap_press & 00269 * 10d0**(7.45d0 * (env_state%temp - const%water_freeze_temp) & 00270 / (env_state%temp - 38d0)) 00271 00272 end function env_state_sat_vapor_pressure 00273 00274 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00275 00276 !> Returns the critical relative humidity from the kappa value (1). 00277 real(kind=dp) function aero_particle_kappa_rh(aero_particle, aero_data, & 00278 env_state) 00279 00280 !> Aerosol particle. 00281 type(aero_particle_t), intent(in) :: aero_particle 00282 !> Aerosol data. 00283 type(aero_data_t), intent(in) :: aero_data 00284 !> Environment state. 00285 type(env_state_t), intent(in) :: env_state 00286 00287 real(kind=dp) :: kappa, diam, C, A 00288 00289 kappa = aero_particle_solute_kappa(aero_particle, aero_data) 00290 A = 4d0 * const%water_surf_eng * const%water_molec_weight & 00291 / (const%univ_gas_const * env_state%temp * const%water_density) 00292 C = sqrt(4d0 * A**3 / 27d0) 00293 diam = vol2diam(aero_particle_volume(aero_particle)) 00294 aero_particle_kappa_rh = C / sqrt(kappa * diam**3) + 1d0 00295 00296 end function aero_particle_kappa_rh 00297 00298 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00299 00300 !> Air density (kg m^{-3}). 00301 real(kind=dp) function env_state_air_den(env_state) 00302 00303 !> Environment state. 00304 type(env_state_t), intent(in) :: env_state 00305 00306 env_state_air_den = const%air_molec_weight & 00307 * env_state_air_molar_den(env_state) 00308 00309 end function env_state_air_den 00310 00311 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00312 00313 !> Air molar density (mol m^{-3}). 00314 real(kind=dp) function env_state_air_molar_den(env_state) 00315 00316 !> Environment state. 00317 type(env_state_t), intent(in) :: env_state 00318 00319 env_state_air_molar_den = env_state%pressure & 00320 / (const%univ_gas_const * env_state%temp) 00321 00322 end function env_state_air_molar_den 00323 00324 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00325 00326 !> Convert (mol m^{-3}) to (ppb). 00327 subroutine gas_state_mole_dens_to_ppb(gas_state, env_state) 00328 00329 !> Gas state. 00330 type(gas_state_t), intent(inout) :: gas_state 00331 !> Environment state. 00332 type(env_state_t), intent(in) :: env_state 00333 00334 gas_state%mix_rat = gas_state%mix_rat & 00335 / env_state_air_molar_den(env_state) * 1d9 00336 00337 end subroutine gas_state_mole_dens_to_ppb 00338 00339 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00340 00341 !> Convert (ppb) to (molecules m^{-3}). 00342 real(kind=dp) function env_state_ppb_to_conc(env_state, ppb) 00343 00344 !> Environment state. 00345 type(env_state_t), intent(in) :: env_state 00346 !> Mixing ratio (ppb). 00347 real(kind=dp), intent(in) :: ppb 00348 00349 env_state_ppb_to_conc = ppb / 1d9 * env_state_air_molar_den(env_state) & 00350 * const%avagadro 00351 00352 end function env_state_ppb_to_conc 00353 00354 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00355 00356 !> Convert (molecules m^{-3}) to (ppb). 00357 real(kind=dp) function env_state_conc_to_ppb(env_state, conc) 00358 00359 !> Environment state. 00360 type(env_state_t), intent(in) :: env_state 00361 !> Concentration (molecules m^{-3}). 00362 real(kind=dp), intent(in) :: conc 00363 00364 env_state_conc_to_ppb = conc * 1d9 / env_state_air_molar_den(env_state) & 00365 / const%avagadro 00366 00367 end function env_state_conc_to_ppb 00368 00369 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00370 00371 !> Do emissions and background dilution from the environment. 00372 subroutine env_state_update_gas_state(env_state, delta_t, & 00373 old_env_state, & 00374 gas_data, gas_state) 00375 00376 !> Current environment. 00377 type(env_state_t), intent(in) :: env_state 00378 !> Time increment to update over. 00379 real(kind=dp), intent(in) :: delta_t 00380 !> Previous environment. 00381 type(env_state_t), intent(in) :: old_env_state 00382 !> Gas data values. 00383 type(gas_data_t), intent(in) :: gas_data 00384 !> Gas state to update. 00385 type(gas_state_t), intent(inout) :: gas_state 00386 00387 real(kind=dp) :: effective_dilution_rate 00388 type(gas_state_t) :: emission, dilution 00389 00390 call gas_state_allocate_size(emission, gas_data%n_spec) 00391 call gas_state_allocate_size(dilution, gas_data%n_spec) 00392 00393 ! account for height changes 00394 effective_dilution_rate = env_state%gas_dilution_rate 00395 if (env_state%height > old_env_state%height) then 00396 effective_dilution_rate = effective_dilution_rate & 00397 + (env_state%height - old_env_state%height) / delta_t / & 00398 old_env_state%height 00399 end if 00400 00401 ! emission = delta_t * gas_emission_rate * gas_emissions 00402 ! but emissions are in (mol m^{-2} s^{-1}) 00403 call gas_state_copy(env_state%gas_emissions, emission) 00404 call gas_state_scale(emission, 1d0 / env_state%height) 00405 call gas_state_mole_dens_to_ppb(emission, env_state) 00406 call gas_state_scale(emission, delta_t * env_state%gas_emission_rate) 00407 00408 ! dilution = delta_t * gas_dilution_rate * (gas_background - gas_state) 00409 call gas_state_copy(env_state%gas_background, dilution) 00410 call gas_state_sub(dilution, gas_state) 00411 call gas_state_scale(dilution, delta_t * effective_dilution_rate) 00412 00413 call gas_state_add(gas_state, emission) 00414 call gas_state_add(gas_state, dilution) 00415 00416 call gas_state_ensure_nonnegative(gas_state) 00417 00418 call gas_state_deallocate(emission) 00419 call gas_state_deallocate(dilution) 00420 00421 end subroutine env_state_update_gas_state 00422 00423 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00424 00425 !> Do emissions and background dilution from the environment for a 00426 !> particle aerosol distribution. 00427 subroutine env_state_update_aero_state(env_state, delta_t, & 00428 old_env_state, bin_grid, aero_data, aero_weight, aero_state, & 00429 n_emit, n_dil_in, n_dil_out) 00430 00431 !> Current environment. 00432 type(env_state_t), intent(in) :: env_state 00433 !> Time increment to update over. 00434 real(kind=dp), intent(in) :: delta_t 00435 !> Previous environment. 00436 type(env_state_t), intent(in) :: old_env_state 00437 !> Bin grid. 00438 type(bin_grid_t), intent(in) :: bin_grid 00439 !> Aero data values. 00440 type(aero_data_t), intent(in) :: aero_data 00441 !> Aero weight. 00442 type(aero_weight_t), intent(in) :: aero_weight 00443 !> Aero state to update. 00444 type(aero_state_t), intent(inout) :: aero_state 00445 !> Number of emitted particles. 00446 integer, intent(out) :: n_emit 00447 !> Number of diluted-in particles. 00448 integer, intent(out) :: n_dil_in 00449 !> Number of diluted-out particles. 00450 integer, intent(out) :: n_dil_out 00451 00452 integer :: i 00453 real(kind=dp) :: sample_prop, effective_dilution_rate 00454 type(aero_state_t) :: aero_state_delta 00455 00456 call aero_state_allocate_size(aero_state_delta, bin_grid%n_bin, & 00457 aero_data%n_spec, aero_data%n_source) 00458 00459 ! account for height changes 00460 effective_dilution_rate = env_state%aero_dilution_rate 00461 if (env_state%height > old_env_state%height) then 00462 effective_dilution_rate = effective_dilution_rate & 00463 + (env_state%height - old_env_state%height) / delta_t / & 00464 old_env_state%height 00465 end if 00466 00467 ! loss to background 00468 sample_prop = 1d0 - exp(- delta_t * effective_dilution_rate) 00469 call aero_state_zero(aero_state_delta) 00470 aero_state_delta%comp_vol = aero_state%comp_vol 00471 call aero_state_sample(aero_state, aero_state_delta, sample_prop, & 00472 AERO_INFO_DILUTION) 00473 n_dil_out = aero_state_total_particles(aero_state_delta) 00474 00475 ! addition from background 00476 sample_prop = 1d0 - exp(- delta_t * effective_dilution_rate) 00477 call aero_state_zero(aero_state_delta) 00478 aero_state_delta%comp_vol = aero_state%comp_vol 00479 call aero_state_add_aero_dist_sample(aero_state_delta, bin_grid, & 00480 aero_data, aero_weight, env_state%aero_background, sample_prop, & 00481 env_state%elapsed_time) 00482 n_dil_in = aero_state_total_particles(aero_state_delta) 00483 call aero_state_add_particles(aero_state, aero_state_delta) 00484 00485 ! emissions 00486 sample_prop = 1d0 & 00487 - exp(- delta_t * env_state%aero_emission_rate / env_state%height) 00488 call aero_state_zero(aero_state_delta) 00489 aero_state_delta%comp_vol = aero_state%comp_vol 00490 call aero_state_add_aero_dist_sample(aero_state_delta, bin_grid, & 00491 aero_data, aero_weight, env_state%aero_emissions, sample_prop, & 00492 env_state%elapsed_time) 00493 n_emit = aero_state_total_particles(aero_state_delta) 00494 call aero_state_add_particles(aero_state, aero_state_delta) 00495 00496 ! update computational volume 00497 aero_state%comp_vol = aero_state%comp_vol * env_state%temp & 00498 / old_env_state%temp 00499 00500 call aero_state_deallocate(aero_state_delta) 00501 00502 end subroutine env_state_update_aero_state 00503 00504 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00505 00506 !> Do emissions and background dilution from the environment for a 00507 !> binned aerosol distribution. 00508 subroutine env_state_update_aero_binned(env_state, delta_t, & 00509 old_env_state, & 00510 bin_grid, aero_data, aero_binned) 00511 00512 !> Current environment. 00513 type(env_state_t), intent(in) :: env_state 00514 !> Time increment to update over. 00515 real(kind=dp), intent(in) :: delta_t 00516 !> Previous environment. 00517 type(env_state_t), intent(in) :: old_env_state 00518 !> Bin grid. 00519 type(bin_grid_t), intent(in) :: bin_grid 00520 !> Aero data values. 00521 type(aero_data_t), intent(in) :: aero_data 00522 !> Aero binned to update. 00523 type(aero_binned_t), intent(inout) :: aero_binned 00524 00525 type(aero_binned_t) :: emission, dilution 00526 real(kind=dp) :: effective_dilution_rate 00527 00528 call aero_binned_allocate_size(emission, bin_grid%n_bin, aero_data%n_spec) 00529 call aero_binned_allocate_size(dilution, bin_grid%n_bin, aero_data%n_spec) 00530 00531 ! account for height changes 00532 effective_dilution_rate = env_state%aero_dilution_rate 00533 if (env_state%height > old_env_state%height) then 00534 effective_dilution_rate = effective_dilution_rate & 00535 + (env_state%height - old_env_state%height) / delta_t / & 00536 old_env_state%height 00537 end if 00538 00539 ! emission = delta_t * aero_emission_rate * aero_emissions 00540 ! but emissions are #/m^2 so we need to divide by height 00541 call aero_binned_add_aero_dist(emission, bin_grid, aero_data, & 00542 env_state%aero_emissions) 00543 call aero_binned_scale(emission, & 00544 delta_t * env_state%aero_emission_rate / env_state%height) 00545 00546 ! dilution = delta_t * aero_dilution_rate 00547 ! * (aero_background - aero_binned) 00548 call aero_binned_add_aero_dist(dilution, bin_grid, aero_data, & 00549 env_state%aero_background) 00550 call aero_binned_sub(dilution, aero_binned) 00551 call aero_binned_scale(dilution, delta_t * effective_dilution_rate) 00552 00553 call aero_binned_add(aero_binned, emission) 00554 call aero_binned_add(aero_binned, dilution) 00555 00556 call aero_binned_deallocate(emission) 00557 call aero_binned_deallocate(dilution) 00558 00559 end subroutine env_state_update_aero_binned 00560 00561 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00562 00563 !> Read environment specification from a spec file. 00564 subroutine spec_file_read_env_state(file, env_state) 00565 00566 !> Spec file. 00567 type(spec_file_t), intent(inout) :: file 00568 !> Environment data. 00569 type(env_state_t), intent(inout) :: env_state 00570 00571 !> \page input_format_env_state Input File Format: Environment State 00572 !! 00573 !! The environment parameters are divided into those specified at 00574 !! the start of the simulation and then either held constant or 00575 !! computed for the rest of the simulation, and those parameters 00576 !! given as prescribed profiles for the entire simulation 00577 !! duration. The variables below are for the first type --- for 00578 !! the prescribed profiles see \ref input_format_env_data. 00579 !! 00580 !! The environment state is specified by the parameters: 00581 !! - \b rel_humidity (real, dimensionless): the relative humidity 00582 !! (0 is completely unsaturated and 1 is fully saturated) 00583 !! - \b pressure (real, unit Pa): the atmospheric pressure 00584 !! - \b latitude (real, unit degrees_north): the latitude of the 00585 !! simulation location 00586 !! - \b longitude (real, unit degrees_east): the longitude of the 00587 !! simulation location 00588 !! - \b altitude (real, unit m): the altitude of the simulation 00589 !! location 00590 !! - \b start_time (real, unit s): the time-of-day of the start of 00591 !! the simulation (in seconds past midnight) 00592 !! - \b start_day (integer): the day-of-year of the start of the 00593 !! simulation (starting from 1 on the first day of the year) 00594 !! 00595 !! See also: 00596 !! - \ref spec_file_format --- the input file text format 00597 !! - \ref output_format_env_state --- the corresponding output 00598 !! format 00599 !! - \ref input_format_env_data --- the prescribed profiles of 00600 !! other environment data 00601 00602 call spec_file_read_real(file, 'rel_humidity', env_state%rel_humid) 00603 call spec_file_read_real(file, 'pressure', env_state%pressure) 00604 call spec_file_read_real(file, 'latitude', env_state%latitude) 00605 call spec_file_read_real(file, 'longitude', env_state%longitude) 00606 call spec_file_read_real(file, 'altitude', env_state%altitude) 00607 call spec_file_read_real(file, 'start_time', env_state%start_time) 00608 call spec_file_read_integer(file, 'start_day', env_state%start_day) 00609 00610 end subroutine spec_file_read_env_state 00611 00612 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00613 00614 !> Average val over all processes. 00615 subroutine env_state_mix(val) 00616 00617 !> Value to average. 00618 type(env_state_t), intent(inout) :: val 00619 00620 #ifdef PMC_USE_MPI 00621 type(env_state_t) :: val_avg 00622 00623 call env_state_allocate(val_avg) 00624 call pmc_mpi_allreduce_average_real(val%temp, val_avg%temp) 00625 call pmc_mpi_allreduce_average_real(val%rel_humid, val_avg%rel_humid) 00626 call pmc_mpi_allreduce_average_real(val%pressure, val_avg%pressure) 00627 val%temp = val_avg%temp 00628 val%rel_humid = val_avg%rel_humid 00629 val%pressure = val_avg%pressure 00630 call env_state_deallocate(val_avg) 00631 #endif 00632 00633 end subroutine env_state_mix 00634 00635 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00636 00637 !> Average val over all processes, with the result only on the root 00638 !> process. 00639 subroutine env_state_reduce_avg(val) 00640 00641 !> Value to average. 00642 type(env_state_t), intent(inout) :: val 00643 00644 #ifdef PMC_USE_MPI 00645 type(env_state_t) :: val_avg 00646 00647 call env_state_allocate(val_avg) 00648 call pmc_mpi_reduce_avg_real(val%temp, val_avg%temp) 00649 call pmc_mpi_reduce_avg_real(val%rel_humid, val_avg%rel_humid) 00650 call pmc_mpi_reduce_avg_real(val%pressure, val_avg%pressure) 00651 if (pmc_mpi_rank() == 0) then 00652 val%temp = val_avg%temp 00653 val%rel_humid = val_avg%rel_humid 00654 val%pressure = val_avg%pressure 00655 end if 00656 call env_state_deallocate(val_avg) 00657 #endif 00658 00659 end subroutine env_state_reduce_avg 00660 00661 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00662 00663 !> Determines the number of bytes required to pack the given value. 00664 integer function pmc_mpi_pack_size_env_state(val) 00665 00666 !> Value to pack. 00667 type(env_state_t), intent(in) :: val 00668 00669 pmc_mpi_pack_size_env_state = & 00670 pmc_mpi_pack_size_real(val%temp) & 00671 + pmc_mpi_pack_size_real(val%rel_humid) & 00672 + pmc_mpi_pack_size_real(val%pressure) & 00673 + pmc_mpi_pack_size_real(val%longitude) & 00674 + pmc_mpi_pack_size_real(val%latitude) & 00675 + pmc_mpi_pack_size_real(val%altitude) & 00676 + pmc_mpi_pack_size_real(val%start_time) & 00677 + pmc_mpi_pack_size_integer(val%start_day) & 00678 + pmc_mpi_pack_size_real(val%elapsed_time) & 00679 + pmc_mpi_pack_size_real(val%solar_zenith_angle) & 00680 + pmc_mpi_pack_size_real(val%height) & 00681 + pmc_mpi_pack_size_gas_state(val%gas_emissions) & 00682 + pmc_mpi_pack_size_real(val%gas_emission_rate) & 00683 + pmc_mpi_pack_size_gas_state(val%gas_background) & 00684 + pmc_mpi_pack_size_real(val%gas_dilution_rate) & 00685 + pmc_mpi_pack_size_aero_dist(val%aero_emissions) & 00686 + pmc_mpi_pack_size_real(val%aero_emission_rate) & 00687 + pmc_mpi_pack_size_aero_dist(val%aero_background) & 00688 + pmc_mpi_pack_size_real(val%aero_dilution_rate) 00689 00690 end function pmc_mpi_pack_size_env_state 00691 00692 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00693 00694 !> Packs the given value into the buffer, advancing position. 00695 subroutine pmc_mpi_pack_env_state(buffer, position, val) 00696 00697 !> Memory buffer. 00698 character, intent(inout) :: buffer(:) 00699 !> Current buffer position. 00700 integer, intent(inout) :: position 00701 !> Value to pack. 00702 type(env_state_t), intent(in) :: val 00703 00704 #ifdef PMC_USE_MPI 00705 integer :: prev_position 00706 00707 prev_position = position 00708 call pmc_mpi_pack_real(buffer, position, val%temp) 00709 call pmc_mpi_pack_real(buffer, position, val%rel_humid) 00710 call pmc_mpi_pack_real(buffer, position, val%pressure) 00711 call pmc_mpi_pack_real(buffer, position, val%longitude) 00712 call pmc_mpi_pack_real(buffer, position, val%latitude) 00713 call pmc_mpi_pack_real(buffer, position, val%altitude) 00714 call pmc_mpi_pack_real(buffer, position, val%start_time) 00715 call pmc_mpi_pack_integer(buffer, position, val%start_day) 00716 call pmc_mpi_pack_real(buffer, position, val%elapsed_time) 00717 call pmc_mpi_pack_real(buffer, position, val%solar_zenith_angle) 00718 call pmc_mpi_pack_real(buffer, position, val%height) 00719 call pmc_mpi_pack_gas_state(buffer, position, val%gas_emissions) 00720 call pmc_mpi_pack_real(buffer, position, val%gas_emission_rate) 00721 call pmc_mpi_pack_gas_state(buffer, position, val%gas_background) 00722 call pmc_mpi_pack_real(buffer, position, val%gas_dilution_rate) 00723 call pmc_mpi_pack_aero_dist(buffer, position, val%aero_emissions) 00724 call pmc_mpi_pack_real(buffer, position, val%aero_emission_rate) 00725 call pmc_mpi_pack_aero_dist(buffer, position, val%aero_background) 00726 call pmc_mpi_pack_real(buffer, position, val%aero_dilution_rate) 00727 call assert(464101191, & 00728 position - prev_position <= pmc_mpi_pack_size_env_state(val)) 00729 #endif 00730 00731 end subroutine pmc_mpi_pack_env_state 00732 00733 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00734 00735 !> Unpacks the given value from the buffer, advancing position. 00736 subroutine pmc_mpi_unpack_env_state(buffer, position, val) 00737 00738 !> Memory buffer. 00739 character, intent(inout) :: buffer(:) 00740 !> Current buffer position. 00741 integer, intent(inout) :: position 00742 !> Value to pack. 00743 type(env_state_t), intent(inout) :: val 00744 00745 #ifdef PMC_USE_MPI 00746 integer :: prev_position 00747 00748 prev_position = position 00749 call pmc_mpi_unpack_real(buffer, position, val%temp) 00750 call pmc_mpi_unpack_real(buffer, position, val%rel_humid) 00751 call pmc_mpi_unpack_real(buffer, position, val%pressure) 00752 call pmc_mpi_unpack_real(buffer, position, val%longitude) 00753 call pmc_mpi_unpack_real(buffer, position, val%latitude) 00754 call pmc_mpi_unpack_real(buffer, position, val%altitude) 00755 call pmc_mpi_unpack_real(buffer, position, val%start_time) 00756 call pmc_mpi_unpack_integer(buffer, position, val%start_day) 00757 call pmc_mpi_unpack_real(buffer, position, val%elapsed_time) 00758 call pmc_mpi_unpack_real(buffer, position, val%solar_zenith_angle) 00759 call pmc_mpi_unpack_real(buffer, position, val%height) 00760 call pmc_mpi_unpack_gas_state(buffer, position, val%gas_emissions) 00761 call pmc_mpi_unpack_real(buffer, position, val%gas_emission_rate) 00762 call pmc_mpi_unpack_gas_state(buffer, position, val%gas_background) 00763 call pmc_mpi_unpack_real(buffer, position, val%gas_dilution_rate) 00764 call pmc_mpi_unpack_aero_dist(buffer, position, val%aero_emissions) 00765 call pmc_mpi_unpack_real(buffer, position, val%aero_emission_rate) 00766 call pmc_mpi_unpack_aero_dist(buffer, position, val%aero_background) 00767 call pmc_mpi_unpack_real(buffer, position, val%aero_dilution_rate) 00768 call assert(205696745, & 00769 position - prev_position <= pmc_mpi_pack_size_env_state(val)) 00770 #endif 00771 00772 end subroutine pmc_mpi_unpack_env_state 00773 00774 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00775 00776 !> Computes the average of val across all processes, storing the 00777 !> result in val_avg on the root process. 00778 subroutine pmc_mpi_reduce_avg_env_state(val, val_avg) 00779 00780 !> Value to average. 00781 type(env_state_t), intent(in) :: val 00782 !> Result. 00783 type(env_state_t), intent(inout) :: val_avg 00784 00785 call env_state_allocate(val_avg) 00786 call env_state_copy(val, val_avg) 00787 call pmc_mpi_reduce_avg_real(val%temp, val_avg%temp) 00788 call pmc_mpi_reduce_avg_real(val%rel_humid, val_avg%rel_humid) 00789 call pmc_mpi_reduce_avg_real(val%pressure, val_avg%pressure) 00790 00791 end subroutine pmc_mpi_reduce_avg_env_state 00792 00793 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00794 00795 !> Write full state. 00796 subroutine env_state_output_netcdf(env_state, ncid) 00797 00798 !> Environment state to write. 00799 type(env_state_t), intent(in) :: env_state 00800 !> NetCDF file ID, in data mode. 00801 integer, intent(in) :: ncid 00802 00803 !> \page output_format_env_state Output File Format: Environment State 00804 !! 00805 !! The environment state NetCDF variables are: 00806 !! - \b temperature (unit K): current air temperature 00807 !! - \b relative_humidity (dimensionless): current air 00808 !! relative humidity (value of 1 means completely saturated) 00809 !! - \b pressure (unit Pa): current air pressure 00810 !! - \b longitude (unit degrees_east): longitude of simulation location 00811 !! - \b latitude (unit degrees_north): latitude of simulation location 00812 !! - \b altitude (unit m): altitude of simulation location 00813 !! - \b start_time_of_day (unit s): time-of-day of the 00814 !! simulation start measured in seconds after midnight UTC 00815 !! - \b start_day_of_year: day-in-year number of the simulation start 00816 !! (starting from 1 on the first day of the year) 00817 !! - \b elapsed_time (unit s): elapsed time since the simulation start 00818 !! - \b solar_zenith_angle (unit radians): current angle from 00819 !! the zenith to the sun 00820 !! - \b height (unit m): current boundary layer mixing height 00821 !! 00822 !! See also: 00823 !! - \ref input_format_env_state and \ref input_format_env_data 00824 !! --- the corresponding input formats 00825 00826 call pmc_nc_write_real(ncid, env_state%temp, "temperature", unit="K", & 00827 standard_name="air_temperature") 00828 call pmc_nc_write_real(ncid, env_state%rel_humid, & 00829 "relative_humidity", unit="1", standard_name="relative_humidity") 00830 call pmc_nc_write_real(ncid, env_state%pressure, "pressure", unit="Pa", & 00831 standard_name="air_pressure") 00832 call pmc_nc_write_real(ncid, env_state%longitude, "longitude", & 00833 unit="degree_east", standard_name="longitude") 00834 call pmc_nc_write_real(ncid, env_state%latitude, "latitude", & 00835 unit="degree_north", standard_name="latitude") 00836 call pmc_nc_write_real(ncid, env_state%altitude, "altitude", unit="m", & 00837 standard_name="altitude") 00838 call pmc_nc_write_real(ncid, env_state%start_time, & 00839 "start_time_of_day", unit="s", description="time-of-day of " & 00840 // "simulation start in seconds since midnight") 00841 call pmc_nc_write_integer(ncid, env_state%start_day, & 00842 "start_day_of_year", & 00843 description="day-of-year number of simulation start") 00844 call pmc_nc_write_real(ncid, env_state%elapsed_time, "elapsed_time", & 00845 unit="s", description="elapsed time since simulation start") 00846 call pmc_nc_write_real(ncid, env_state%solar_zenith_angle, & 00847 "solar_zenith_angle", unit="radian", & 00848 description="current angle from the zenith to the sun") 00849 call pmc_nc_write_real(ncid, env_state%height, "height", unit="m", & 00850 long_name="boundary layer mixing height") 00851 00852 end subroutine env_state_output_netcdf 00853 00854 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00855 00856 !> Read full state. 00857 subroutine env_state_input_netcdf(env_state, ncid) 00858 00859 !> Environment state to read. 00860 type(env_state_t), intent(inout) :: env_state 00861 !> NetCDF file ID, in data mode. 00862 integer, intent(in) :: ncid 00863 00864 call pmc_nc_read_real(ncid, env_state%temp, "temperature") 00865 call pmc_nc_read_real(ncid, env_state%rel_humid, "relative_humidity") 00866 call pmc_nc_read_real(ncid, env_state%pressure, "pressure") 00867 call pmc_nc_read_real(ncid, env_state%longitude, "longitude") 00868 call pmc_nc_read_real(ncid, env_state%latitude, "latitude") 00869 call pmc_nc_read_real(ncid, env_state%altitude, "altitude") 00870 call pmc_nc_read_real(ncid, env_state%start_time, & 00871 "start_time_of_day") 00872 call pmc_nc_read_integer(ncid, env_state%start_day, & 00873 "start_day_of_year") 00874 call pmc_nc_read_real(ncid, env_state%elapsed_time, "elapsed_time") 00875 call pmc_nc_read_real(ncid, env_state%solar_zenith_angle, & 00876 "solar_zenith_angle") 00877 call pmc_nc_read_real(ncid, env_state%height, "height") 00878 00879 end subroutine env_state_input_netcdf 00880 00881 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00882 00883 end module pmc_env_state