def water_viscosity_kinematic(temperature=25*unit('degC'),pressure=1*unit('atm')): ''' Return the kinematic viscosity of water in m2/s = Stokes at the specified temperature. Parameters ---------- temperature : Quantity, optional The temperature. Defaults to 25 degC if omitted. pressure : Quantity, optional The ambient pressure of the solution. Defaults to atmospheric pressure (1 atm) if omitted. Returns ------- Quantity The kinematic viscosity of water in Stokes (m2/s) Examples -------- >>> water_viscosity_kinematic() #doctest: +ELLIPSIS <Quantity(8.899146003595295e-07, 'meter ** 2 / second')> See Also -------- water_viscosity_dynamic water_density ''' kviscosity = water_viscosity_dynamic(temperature,pressure) / water_density(temperature,pressure) logger.info('Computed kinematic viscosity of water as %s at T=%s and P = %s ' % (kviscosity,temperature,pressure)) return kviscosity.to('m**2 / s')
def water_specific_weight(temperature=25*unit('degC'),pressure=1*unit('atm')): ''' Return the specific weight of water in N/m3 at the specified temperature and pressure. Parameters ---------- temperature : Quantity, optional The temperature. Defaults to 25 degC if omitted. pressure : Quantity, optional The ambient pressure of the solution. Defaults to atmospheric pressure (1 atm) if omitted. Returns ------- Quantity The specific weight of water in N/m3. Examples -------- >>> water_specific_weight() #doctest: +ELLIPSIS <Quantity(9777.637025975, 'newton / meter ** 3')> See Also -------- water_density ''' spweight = water_density(temperature,pressure) * unit.g_n logger.info('Computed specific weight of water as %s at T=%s and P = %s' % (spweight,temperature,pressure)) return spweight.to('N/m ** 3')
def _debye_parameter_activity(temperature="25 degC"): """ Return the constant A for use in the Debye-Huckel limiting law (base 10) Parameters ---------- temperature : str Quantity, optional String representing the temperature of the solution. Defaults to '25 degC' if not specified. Returns ------- Quantity The parameter A for use in the Debye-Huckel limiting law (base e) Notes ----- The parameter A is equal to: [#]_ .. math:: A^{\\gamma} = {e^3 ( 2 \\pi N_A {\\rho})^{0.5} \\over (4 \\pi \\epsilon_o \\epsilon_r k T)^{1.5}} Note that this equation returns the parameter value that can be used to calculate the natural logarithm of the activity coefficient. For base 10, divide the value returned by 2.303. The value is often given in base 10 terms (0.509 at 25 degC) in older textbooks. References ---------- .. [#] Archer, Donald G. and Wang, Peiming. "The Dielectric Constant of Water \ and Debye-Huckel Limiting Law Slopes." /J. Phys. Chem. Ref. Data/ 19(2), 1990. Examples -------- >>> _debye_parameter_activity() #doctest: +ELLIPSIS 1.17499... See Also -------- _debye_parameter_osmotic """ debyeparam = ( unit.elementary_charge ** 3 * (2 * math.pi * unit.avogadro_number * h2o.water_density(unit(temperature))) ** 0.5 / ( 4 * math.pi * unit.epsilon_0 * h2o.water_dielectric_constant(unit(temperature)) * unit.boltzmann_constant * unit(temperature) ) ** 1.5 ) logger.info("Computed Debye-Huckel Limiting Law Constant A^{\\gamma} = %s at %s" % (debyeparam, temperature)) return debyeparam.to("kg ** 0.5 / mol ** 0.5")
def water_dielectric_constant(temperature=25*unit('degC')): ''' Return the dielectric constant of water at the specified temperature. Parameters ---------- temperature : Quantity, optional The temperature. Defaults to 25 degC if omitted. Returns ------- float The dielectric constant (or permittivity) of water relative to the permittivity of a vacuum. Dimensionless. Notes ----- This function implements a quadratic fit of measured permittivity data as reported in the CRC Handbook [#]_. The parameters given are valid over the range 273 K to 372 K. Permittivity should not be extrapolated beyond this range. .. math:: \\epsilon(T) = a + b T + c T^2 References ---------- .. [#] "Permittivity (Dielectric Constant) of Liquids." *CRC Handbook of Chemistry and Physics*, 92nd ed, pp 6-187 - 6-208. Examples -------- >>> water_dielectric_constant(unit('20 degC')) #doctest: +ELLIPSIS 80.15060... Display an error if 'temperature' is outside the valid range >>> water_dielectric_constant(-5*unit('degC')) ''' # do not return anything if 'temperature' is outside the range for which # this fit applies if temperature < 273 * unit('K') or temperature > 372 * unit('K'): logger.error('Specified temperature (%s) exceeds valid range of data. Cannot extrapolate.' % temperature.to('K')) return None # otherwise, calculate the dielectric constant using the quadratic fit a = 0.24921e3 b = -0.79069e0 c = 0.72997e-3 dielectric = a + b * temperature.to('K').magnitude + c * temperature.to('K').magnitude ** 2 logger.info('Computed dielectric constant of water as %s at %s' % (dielectric,temperature)) logger.debug('Computed dielectric constant of water using empirical equation given in "Permittivity (Dielectric Constant) of Liquids." CRC Handbook of Chemistry and Physics, 92nd ed, pp 6-187 - 6-208.') return dielectric
def test_molar_conductivity_water(self): self.s1 = pyEQL.Solution(temperature='25 degC') actual = self.s1.get_molar_conductivity('H2O').to('m**2*S/mol') expected = pyEQL.unit('0 m**2 * S / mol') result = (actual-expected).magnitude self.assertAlmostEqual(result,0,self.places)
def test_molar_conductivity_neutral(self): self.s1 = pyEQL.Solution([['FeCl3','0.001 mol/L']],temperature='25 degC') actual = self.s1.get_molar_conductivity('FeCl3').to('m**2*S/mol') expected = pyEQL.unit('0 m**2 * S / mol') result = (actual-expected).magnitude self.assertAlmostEqual(result,0,self.places)
def test_molar_conductivity_chloride(self): # Cl- - 76.31 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution([['Na+','0.001 mol/L'],['Cl-','0.001 mol/L']],temperature='25 degC') result = self.s1.get_molar_conductivity('Cl-').to('m**2*S/mol').magnitude expected = pyEQL.unit('76.31e-4 m**2 * S / mol').magnitude self.assertWithinExperimentalError(result,expected,self.tol)
def test_molar_conductivity_magnesium(self): # Mg+2 - 106 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution([['Mg+2','0.001 mol/L'],['Cl-','0.002 mol/L']],temperature='25 degC') result = self.s1.get_molar_conductivity('Mg+2').to('m**2*S/mol').magnitude expected = pyEQL.unit('106e-4 m**2 * S / mol').magnitude self.assertWithinExperimentalError(result,expected,self.tol)
def test_molar_conductivity_potassium(self): # K+ - 73.48 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution([['K+','0.001 mol/L'],['Cl-','0.001 mol/L']],temperature='25 degC') result = self.s1.get_molar_conductivity('K+').to('m**2*S/mol').magnitude expected = pyEQL.unit('73.48e-4 m**2 * S / mol').magnitude self.assertWithinExperimentalError(result,expected,self.tol)
def test_molar_conductivity_hydrogen(self): # H+ - 349.65 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution(temperature='25 degC') result = self.s1.get_molar_conductivity('H+').to('m**2*S/mol').magnitude expected = pyEQL.unit('349.65e-4 m**2 * S / mol').magnitude self.assertWithinExperimentalError(result,expected,self.tol)
def test_molar_conductivity_hydroxide(self): # OH- - 198 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution(temperature='25 degC') result = self.s1.get_molar_conductivity('OH-').to('m**2*S/mol').magnitude expected = pyEQL.unit('198e-4 m**2 * S / mol').magnitude self.assertWithinExperimentalError(result,expected,self.tol)
def test_molar_conductivity_sulfate(self): # SO4-2 - 160 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution([['Na+','0.002 mol/L'],['SO4-2','0.001 mol/L']],temperature='25 degC') result = self.s1.get_molar_conductivity('SO4-2').to('m**2*S/mol').magnitude expected = pyEQL.unit('160.0e-4 m**2 * S / mol').magnitude self.assertWithinExperimentalError(result,expected,self.tol)
def test_effective_pitzer_mgcl2_activity(self): # test the activity coefficient of MgCl2 # corresponds to 0.515m, 1.03m, 2.58m, and 4.1m multiple = [1,2,5,8] expected=[0.5,0.5,0.67,1.15] # import the parameters database from pyEQL import paramsDB as db for item in range(len(multiple)): s1 = self.mock_seawater(multiple[item]) Salt = pyEQL.salt_ion_match.Salt('Mg+2','Cl-') db.search_parameters(Salt.formula) param = db.get_parameter(Salt.formula,'pitzer_parameters_activity') alpha1 = 2 alpha2 = 0 molality = Salt.get_effective_molality(s1.get_ionic_strength()) temperature = str(s1.get_temperature()) activity_coefficient=pyEQL.activity_correction.get_activity_coefficient_pitzer(s1.get_ionic_strength(), \ molality,alpha1,alpha2,param.get_value()[0],param.get_value()[1],param.get_value()[2],param.get_value()[3], \ Salt.z_cation,Salt.z_anion,Salt.nu_cation,Salt.nu_anion,temperature) # convert the result to a rational activity coefficient result = activity_coefficient * (1+pyEQL.unit('0.018 kg/mol')*s1.get_total_moles_solute()/s1.get_solvent_mass()) #print(result,expected[item]) self.assertWithinExperimentalError(result,expected[item],self.tol)
def _debye_parameter_B(temperature='25 degC'): ''' Return the constant B used in the extended Debye-Huckel equation Parameters ---------- temperature : str Quantity, optional String representing the temperature of the solution. Defaults to '25 degC' if not specified. Notes ----- The parameter B is equal to: [#]_ .. math:: B = ( {8 \\pi N_A e^2 \\over 1000 \\epsilon k T} ) ^ {1 \\over 2} .. [#] Bockris and Reddy. /Modern Electrochemistry/, vol 1. Plenum/Rosetta, 1977, p.210. Examples -------- >>> _debye_parameter_B() #doctest: +ELLIPSIS 0.3291... ''' # TODO - fix this and resolve units param_B = ( 8 * math.pi * unit.avogadro_number * unit.elementary_charge ** 2 / (h2o.water_density(unit(temperature)) * unit.epsilon_0 * h2o.water_dielectric_constant(unit(temperature)) * unit.boltzmann_constant * unit(temperature)) )** 0.5 return param_B.to_base_units()
def test_molar_conductivity_hydrogen(self): # H+ - 349.65 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution(temperature='25 degC') actual = self.s1.get_molar_conductivity('H+').to('m**2*S/mol') expected = pyEQL.unit('349.65e-4 m**2 * S / mol') result = (actual-expected).magnitude self.assertAlmostEqual(result,0,self.places)
def test_molar_conductivity_sulfate(self): # SO4-2 - 160 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution([['Na+','0.002 mol/L'],['SO4-2','0.001 mol/L']],temperature='25 degC') actual = self.s1.get_molar_conductivity('SO4-2').to('m**2*S/mol') expected = pyEQL.unit('160.0e-4 m**2 * S / mol') result = (actual-expected).magnitude self.assertAlmostEqual(result,0,self.places)
def test_molar_conductivity_magnesium(self): # Mg+2 - 106 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution([['Mg+2','0.001 mol/L'],['Cl-','0.002 mol/L']],temperature='25 degC') actual = self.s1.get_molar_conductivity('Mg+2').to('m**2*S/mol') expected = pyEQL.unit('106e-4 m**2 * S / mol') result = (actual-expected).magnitude self.assertAlmostEqual(result,0,self.places)
def test_molar_conductivity_chloride(self): # Cl- - 76.31 x 10 ** -4 m ** 2 S / mol self.s1 = pyEQL.Solution([['Na+','0.001 mol/L'],['Cl-','0.001 mol/L']],temperature='25 degC') actual = self.s1.get_molar_conductivity('Cl-').to('m**2*S/mol') expected = pyEQL.unit('76.31e-4 m**2 * S / mol') result = (actual-expected).magnitude self.assertAlmostEqual(result,0,self.places)
def water_density(temperature=25*unit('degC'),pressure=1*unit('atm')): # TODO add pressure?? # TODO more up to date equation?? ''' Return the density of water in kg/m3 at the specified temperature and pressure. Parameters ---------- temperature : float or int, optional The temperature in Celsius. Defaults to 25 degrees if not specified. pressure : float or int, optional The ambient pressure of the solution in Pascals (N/m2). Defaults to atmospheric pressure (101325 Pa) if not specified. Returns ------- float The density of water in kg/m3. Notes ----- Based on the following empirical equation reported in [#]_ .. math:: \\rho_W = 999.65 + 0.20438 T - 6.1744e-2 T ^ {1.5} Where :math:`T` is the temperature in Celsius. .. [#] Sohnel, O and Novotny, P. *Densities of Aqueous Solutions of Inorganic Substances.* Elsevier Science, Amsterdam, 1985. Examples -------- >>> water_density(25*unit('degC')) #doctest: +ELLIPSIS <Quantity(997.0415, 'kilogram / meter ** 3')> ''' # calculate the magnitude density = 999.65 + 0.20438 * temperature.to('degC').magnitude - 6.1744e-2 * temperature.to('degC').magnitude ** 1.5 # assign the proper units density = density * unit('kg/m**3') logger.info('Computed density of water as %s at T= %s and P = %s' % (density,temperature,pressure)) logger.debug('Computed density of water using empirical relation in Sohnel and Novotny, "Densities of Aqueous Solutions of Inorganic Substances," 1985' ) return density.to('kg/m**3')
def __init__(self,formula,amount,volume,solvent_mass,parameters={}): ''' Parameters ---------- formula : str Chemical formula for the solute. Charged species must contain a + or - and (for polyvalent solutes) a number representing the net charge (e.g. 'SO4-2'). amount : str The amount of substance in the specified unit system. The string should contain both a quantity and a pint-compatible representation of a unit. e.g. '5 mol/kg' or '0.1 g/L' volume : pint Quantity The volume of the solution solvent_mass : pint Quantity The mass of solvent in the parent solution. parameters : dictionary, optional Dictionary of custom parameters, such as diffusion coefficients, transport numbers, etc. Specify parameters as key:value pairs separated by commas within curly braces, e.g. {diffusion_coeff:5e-10,transport_number:0.8}. The 'key' is the name that will be used to access the parameter, the value is its value. ''' # import the chemical formula interpreter module import pyEQL.chemical_formula as chem # check that 'formula' is a valid chemical formula if not chem.is_valid_formula: logger.error('Invalid chemical formula specified.') return None else: self.formula = formula # set molecular weight self.mw = chem.get_molecular_weight(formula) * unit('g/mol') # set formal charge self.charge = chem.get_formal_charge(formula) # translate the 'amount' string into a pint Quantity quantity = unit(amount) self.moles = quantity.to('moles','chem',mw=self.mw,volume=volume,solvent_mass=solvent_mass) # trigger the function that checks whether parameters already exist for this species, and if not, # searches the database files and creates them db.search_parameters(self.formula)
def _debye_parameter_volume(temperature='25 degC'): ''' Return the constant A_V, the Debye-Huckel limiting slope for apparent molar volume. Parameters ---------- temperature : str Quantity, optional String representing the temperature of the solution. Defaults to '25 degC' if not specified. Notes ----- Takes the value 1.8305 cm ** 3 * kg ** 0.5 / mol ** 1.5 at 25 C. This constant is calculated according to: [#]_ .. math:: A_V = -2 A_{\\phi} R T [ {3 \\over \\epsilon} {{\\partial \\epsilon \\over \\partial p} \ } - {{1 \\over \\rho}{\\partial \\rho \\over \\partial p} }] NOTE: at this time, the term in brackets (containing the partial derivatives) is approximate. These approximations give the correct value of the slope at 25 degC and produce estimates with less than 10% error between 0 and 60 degC. The derivative of epsilon with respect to pressure is assumed constant (for atmospheric pressure) at -0.01275 1/MPa. Note that the negative sign does not make sense in light of real data, but is required to give the correct result. The second term is equivalent to the inverse of the bulk modulus of water, which is taken to be 2.2 GPa. [#]_ References ---------- .. [#] Archer, Donald G. and Wang, Peiming. "The Dielectric Constant of Water \ and Debye-Huckel Limiting Law Slopes." /J. Phys. Chem. Ref. Data/ 19(2), 1990. .. [#] http://hyperphysics.phy-astr.gsu.edu/hbase/permot3.html Examples -------- TODO See Also -------- _debye_parameter_osmotic ''' # TODO - add partial derivatives to calculation epsilon = h2o.water_dielectric_constant(unit(temperature)) dedp = unit('-0.01275 1/MPa') result = -2 * _debye_parameter_osmotic(temperature) * unit.R * unit(temperature) * \ (3 / epsilon * dedp - 1/unit('2.2 GPa')) #result = unit('1.898 cm ** 3 * kg ** 0.5 / mol ** 1.5') if unit(temperature) != unit('25 degC'): logger.warning('Debye-Huckel limiting slope for volume is approximate when T is not equal to 25 degC') logger.info('Computed Debye-Huckel Limiting Slope for volume A^V = %s at %s' % (result,temperature)) return result.to('cm ** 3 * kg ** 0.5 / mol ** 1.5')
def set_moles(self,amount,volume,solvent_mass): ''' Set the amount of a substance present in the solution Parameters ---------- amount: str quantity Desired amount of substance. Must be greater than or equal to zero and given in mass or substance units. ''' quantity = unit(amount) self.moles = quantity.to('moles','chem',mw=self.mw,volume=volume,solvent_mass=solvent_mass)
def add_moles(self,amount,volume,solvent_mass): ''' Increase or decrease the amount of a substance present in the solution Parameters ---------- amount: str quantity Amount of substance to add. Must be in mass or substance units. Negative values indicate subtraction of material. ''' quantity = unit(amount) self.moles += quantity.to('moles','chem',mw=self.mw,volume=volume,solvent_mass=solvent_mass)
def get_activity_coefficient_davies(ionic_strength, formal_charge=1, temperature="25 degC"): """ Return the activity coefficient of solute in the parent solution according to the Davies equation. Parameters ---------- formal_charge : int, optional The charge on the solute, including sign. Defaults to +1 if not specified. ionic_strength : Quantity The ionic strength of the parent solution, mol/kg temperature : str Quantity, optional String representing the temperature of the solution. Defaults to '25 degC' if not specified. Returns ------- Quantity The mean molal (mol/kg) scale ionic activity coefficient of solute, dimensionless. See Also -------- _debye_parameter_activity Notes ----- Activity coefficient is calculated according to: [#]_ .. math:: \\ln \\gamma = A^{\\gamma} z_i^2 ({\sqrt I \\over (1 + \sqrt I)} + 0.2 I) Valid for 0.1 < I < 0.5 References ---------- .. [#] Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 103. Wiley Interscience, 1996. """ # check if this method is valid for the given ionic strength if not ionic_strength.magnitude <= 0.5 and ionic_strength.magnitude >= 0.1: logger.warning("Ionic strength exceeds valid range of the Davies equation") # the units in this empirical equation don't work out, so we must use magnitudes log_f = ( -_debye_parameter_activity(temperature).magnitude * formal_charge ** 2 * (ionic_strength.magnitude ** 0.5 / (1 + ionic_strength.magnitude ** 0.5) - 0.2 * ionic_strength.magnitude) ) return math.exp(log_f) * unit("1 dimensionless")
def get_activity_coefficient_debyehuckel(ionic_strength, formal_charge=1, temperature="25 degC"): """ Return the activity coefficient of solute in the parent solution according to the Debye-Huckel limiting law. Parameters ---------- formal_charge : int, optional The charge on the solute, including sign. Defaults to +1 if not specified. ionic_strength : Quantity The ionic strength of the parent solution, mol/kg temperature : str Quantity, optional String representing the temperature of the solution. Defaults to '25 degC' if not specified. Returns ------- Quantity The mean molal (mol/kg) scale ionic activity coefficient of solute, dimensionless. See Also -------- _debye_parameter_activity Notes ----- Activity coefficient is calculated according to: [#]_ .. math:: \\ln \\gamma = A^{\\gamma} z_i^2 \sqrt I Valid only for I < 0.005 References ---------- .. [#] Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 103. Wiley Interscience, 1996. """ # check if this method is valid for the given ionic strength if not ionic_strength.magnitude <= 0.005: logger.warning("Ionic strength exceeds valid range of the Debye-Huckel limiting law") log_f = -_debye_parameter_activity(temperature) * formal_charge ** 2 * ionic_strength ** 0.5 return math.exp(log_f) * unit("1 dimensionless")
def get_activity_coefficient_guntelberg(ionic_strength,formal_charge=1,temperature='25 degC'): ''' Return the activity coefficient of solute in the parent solution according to the Guntelberg approximation. Parameters ---------- formal_charge : int, optional The charge on the solute, including sign. Defaults to +1 if not specified. ionic_strength : Quantity The ionic strength of the parent solution, mol/kg temperature : str Quantity, optional String representing the temperature of the solution. Defaults to '25 degC' if not specified. Returns ------- Quantity The mean molal (mol/kg) scale ionic activity coefficient of solute, dimensionless. See Also -------- _debye_parameter_activity Notes ------ Activity coefficient is calculated according to: [#]_ .. math:: \\ln \\gamma = A^{\\gamma} z_i^2 {\sqrt I \\over (1 + \sqrt I)} Valid for I < 0.1 References ---------- .. [#] Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 103. Wiley Interscience, 1996. ''' # check if this method is valid for the given ionic strength if not ionic_strength.magnitude <= 0.1: logger.warning('Ionic strength exceeds valid range of the Guntelberg approximation') log_f = - _debye_parameter_activity(temperature) * formal_charge ** 2 * ionic_strength ** 0.5 / (1+ionic_strength.magnitude ** 0.5) return math.exp(log_f) * unit('1 dimensionless')
def adjust_temp_pitzer(c1,c2,c3,c4,c5,temp,temp_ref=unit('298.15 K')): ''' Calculate a parameter for th e Pitzer model based on temperature-dependent coefficients c1,c2,c3,c4,and c5. Parameters ---------- c1, c2, c3, c4, c5: float Temperature-dependent coefficients for the pitzer parameter of interest. temp: Quantity The temperature at which the Pitzer parameter is to be calculated temp_ref: Quantity, optional The reference temperature on which the parameters are based. 298.15 K if omitted. As described in the PHREEQC documentation ''' pitzer_param = c1 + c2 * (1/temp + 1/temp_ref) + c2 * math.log(temp/temp_ref) \ + c3 * (temp - temp_ref) + c4 * (temp ** 2 - temp_ref ** 2) + c5 * (temp ** -2 - temp_ref ** -2) return pitzer_param
def _debye_parameter_volume(temperature='25 degC'): ''' Return the constant A_V, the Debye-Huckel limiting slope for apparent molar volume. Parameters ---------- temperature : str Quantity, optional String representing the temperature of the solution. Defaults to '25 degC' if not specified. Notes ----- Takes the value 1.8305 cm ** 3 * kg ** 0.5 / mol ** 1.5 at 25 C. This constant is calculated according to: [#]_ .. math:: A_V = -2 A_{\\phi} R T [ {3 \\over \\epsilon} {{\\partial \\epsilon \\over \\partial p} \ } - {{1 \\over \\rho}{\\partial \\rho \\over \\partial p} }] NOTE: at this time, the term in brackets (containing the partial derivatives) is approximate. These approximations give the correct value of the slope at 25 degC and produce estimates with less than 10% error between 0 and 60 degC. The derivative of epsilon with respect to pressure is assumed constant (for atmospheric pressure) at -0.01275 1/MPa. Note that the negative sign does not make sense in light of real data, but is required to give the correct result. The second term is equivalent to the inverse of the bulk modulus of water, which is taken to be 2.2 GPa. [#]_ References ---------- .. [#] Archer, Donald G. and Wang, Peiming. "The Dielectric Constant of Water \ and Debye-Huckel Limiting Law Slopes." /J. Phys. Chem. Ref. Data/ 19(2), 1990. .. [#] http://hyperphysics.phy-astr.gsu.edu/hbase/permot3.html Examples -------- TODO See Also -------- _debye_parameter_osmotic ''' # TODO - add partial derivatives to calculation epsilon = h2o.water_dielectric_constant(unit(temperature)) dedp = unit('-0.01275 1/MPa') result = -2 * _debye_parameter_osmotic(temperature) * unit.R * unit(temperature) * \ (3 / epsilon * dedp - 1/unit('2.2 GPa')) #result = unit('1.898 cm ** 3 * kg ** 0.5 / mol ** 1.5') if unit(temperature) != unit('25 degC'): logger.warning( 'Debye-Huckel limiting slope for volume is approximate when T is not equal to 25 degC' ) logger.info( 'Computed Debye-Huckel Limiting Slope for volume A^V = %s at %s' % (result, temperature)) return result.to('cm ** 3 * kg ** 0.5 / mol ** 1.5')
def water_viscosity_dynamic(temperature=25*unit('degC'),pressure=1*unit('atm')): ''' Return the dynamic (absolute) viscosity of water in N-s/m2 = Pa-s = kg/m-s at the specified temperature. Parameters ---------- temperature : Quantity, optional The temperature. Defaults to 25 degC if omitted. pressure : Quantity, optional The ambient pressure of the solution. Defaults to atmospheric pressure (1 atm) if omitted. Returns ------- Quantity The dynamic (absolute) viscosity of water in N-s/m2 = Pa-s = kg/m-s Notes ----- Implements the international equation for viscosity of water as specified by NIST [#]_ Valid for 273 < temperature < 1073 K and 0 < pressure < 100,000,000 Pa References ---------- .. [#] Sengers, J.V. "Representative Equations for the Viscosity of Water Substance." *J. Phys. Chem. Ref. Data* 13(1), 1984.http://www.nist.gov/data/PDFfiles/jpcrd243.pdf Examples -------- >>> water_viscosity_dynamic(20*unit('degC')) #doctest: +ELLIPSIS <Quantity(0.000998588610804179, 'kilogram / meter / second')> >>> water_viscosity_dynamic(unit('100 degC'),unit('25 MPa')) #doctest: +ELLIPSIS <Quantity(0.00028165034364318573, 'kilogram / meter / second')> >>> water_viscosity_dynamic(25*unit('degC'),0.1*unit('MPa')) #doctest: +ELLIPSIS <Quantity(0.0008872817880143659, 'kilogram / meter / second')> #TODO - check these again after I implement pressure-dependent density function ''' # generate warnings if temp or pressure are outside valid range of equation if temperature < 273 * unit('K') or temperature > 1073 * unit('K'): logger.error('Specified temperature (%s) exceeds valid range of NIST equation for viscosity of water. Cannot extrapolate.' % temperature) return None if pressure < 0 * unit('Pa') or pressure > 100000000 * unit ('Pa'): logger.error('Specified pressure (%s) exceeds valid range of NIST equation for viscosity of water. Cannot extrapolate.' % pressure) return None # calculate dimensionless temperature and pressure T_star = 647.27 #K P_star = 22115000 #Pa rho_star = 317.763 #kg/m3 T_bar = temperature.to('K').magnitude / T_star P_bar = pressure.to('Pa').magnitude / P_star rho_bar = water_density(temperature,pressure).magnitude / rho_star # calculate the first function, mu_o mu_star = 1e-6 #Pa-s a = [0.0181583,0.0177624,0.0105287,-0.0036477] sum_o = 0 mu_temp = 0 for index in range(len(a)): sum_o += a[index] * T_bar ** -index mu_o = mu_star * math.sqrt(T_bar) / sum_o # calculate the second fucntion, mu_1 b=[[0.501938,0.235622,-0.274637,0.145831,-0.0270448],[0.162888,0.789393,-0.743539,0.263129,-0.0253093],[-0.130356,0.673665,-0.959456,0.347247,-0.0267758],[0.907919,1.207552,-0.687343,0.213486,-0.0822904],[-0.551119,0.0670665,-0.497089,0.100754,0.0602253],[0.146543,-0.0843370,0.195286,-0.032932,-0.0202595]] mu_1 = 0 for i in range(len(b)): for j in range(len(b[i])): mu_temp += rho_bar * b[i][j] * (1/T_bar -1 ) ** i * (rho_bar -1) ** j mu_1 = math.exp(mu_temp) # multiply the functions to return the viscosity viscosity = mu_o * mu_1 *unit('kg/m/s') logger.info('Computed dynamic (absolute) viscosity of water as %s at T=%s and P = %s' % (viscosity,temperature,pressure)) logger.debug('Computed dynamic (absolute) viscosity of water using empirical NIST equation described in Sengers, J.V. "Representative Equations for the Viscosity of Water Substance." J. Phys. Chem. Ref. Data 13(1), 1984.') return viscosity.to('kg/m/s')
def donnan_eql(solution,fixed_charge): ''' Return a solution object in equilibrium with fixed_charge Parameters ---------- Solution : Solution object The external solution to be brought into equilibrium with the fixed charges fixed_charge : str quantity String representing the concentration of fixed charges, including sign. May be specified in mol/L or mol/kg units. e.g. '1 mol/kg' Returns ------- Solution A solution that has established Donnan equilibrium with the external (input) Solution Notes ----- The general equation representing the equilibrium between an external electrolyte solution and an ion-exchange medium containing fixed charges is:[#]_ .. math:: {a_- \\over \\bar a_-}^{1 \\over z_-} {\\bar a_+ \\over a_+}^{1 \\over z_+} \ = exp({\\Delta \\pi \\bar V \\over {RT z_+ \\nu_+}}) Where subscripts :math:`+` and :math:`-` indicate the cation and anion, respectively, the overbar indicates the membrane phase, :math:`a` represents activity, :math:`z` represents charge, :math:`\\nu` represents the stoichiometric coefficient, :math:`V` represents the partial molar volume of the salt, and :math:`\\Delta \\pi` is the difference in osmotic pressure between the membrane and the solution phase. In addition, electroneutrality must prevail within the membrane phase: .. math:: \\bar C_+ z_+ + \\bar X + \\bar C_- z_- = 0 Where :math:`C` represents concentration and :math:`X` is the fixed charge concentration in the membrane or ion exchange phase. This function solves these two equations simultaneously to arrive at the concentrations of the cation and anion in the membrane phase. It returns a solution equal to the input solution except that the concentrations of the predominant cation and anion have been adjusted according to this equilibrium. NOTE that this treatment is only capable of equilibrating a single salt. This salt is identified by the get_salt() method. References ---------- .. [#] Strathmann, Heiner, ed. *Membrane Science and Technology* vol. 9, 2004. \ Chapter 2, p. 51. http://dx.doi.org/10.1016/S0927-5193(04)80033-0 Examples -------- TODO See Also -------- get_salt() ''' # identify the salt salt = solution.get_salt() # convert fixed_charge in to a quantity fixed_charge = unit(fixed_charge) # identify variables from the external solution conc_cation_soln = solution.get_amount(salt.cation,str(fixed_charge.units)) conc_anion_soln = solution.get_amount(salt.anion,str(fixed_charge.units)) act_cation_soln = solution.get_activity(salt.cation) act_anion_soln = solution.get_activity(salt.anion) z_cation= salt.z_cation z_anion = salt.z_anion nu_cation = salt.nu_cation # get the partial molar volume for the salt, or calculate it from the ions # TODO - consider how to incorporate pitzer parameters if db.has_parameter(salt.formula,'partial_molar_volume'): item = db.get_parameter(salt.formula,'partial_molar_volume') molar_volume = item.get_value() elif db.has_parameter(salt.cation,'partial_molar_volume') and db.has_parameter(salt.anion,'partial_molar_volume'): cation_vol = solution.get_solute(salt.cation).get_parameter('partial_molar_volume') anion_vol = solution.get_solute(salt.anion).get_parameter('partial_molar_volume') molar_volume = cation_vol + anion_vol else: logger.error('Required partial molar volume information not available. Aborting.') return None # initialize the equilibrated solution - start with a direct copy of the # input / external solution donnan_soln = solution.copy() # do nothing if either of the ion concentrations is zero if conc_cation_soln.magnitude == 0 or conc_anion_soln.magnitude == 0: return donnan_soln # define a function representing the donnan equilibrium as a function # of the two unknown actvities to feed to the nonlinear solver # the stuff in the term below doesn't change on iteration, so calculate it up-front # assign it the correct units and extract the magnitude for a performance gain exp_term = (molar_volume / (unit.R * solution.get_temperature() * z_cation * nu_cation)).to('1/Pa').magnitude def donnan_solve(x): '''Where x is the magnitude of co-ion concentration ''' # solve for the counter-ion concentration by enforcing electroneutrality # using only floats / ints here instead of quantities helps performance if fixed_charge.magnitude >= 0: # counter-ion is the anion conc_cation_mem = x / abs(z_cation) conc_anion_mem = -(conc_cation_mem * z_cation + fixed_charge.magnitude) / z_anion elif fixed_charge.magnitude < 0: # counter-ion is the cation conc_anion_mem = x / abs(z_anion) conc_cation_mem = -(conc_anion_mem * z_anion + fixed_charge.magnitude) / z_cation # match the units given for fixed_charge units = str(fixed_charge.units) # set the cation and anion concentrations in the membrane phase equal # to the current guess donnan_soln.set_amount(salt.cation,str(conc_cation_mem)+units) donnan_soln.set_amount(salt.anion,str(conc_anion_mem)+units) # get the new concentrations and activities act_cation_mem = donnan_soln.get_activity(salt.cation) act_anion_mem = donnan_soln.get_activity(salt.anion) # compute the difference in osmotic pressure # using the magnitudes here helps performance delta_pi = donnan_soln.get_osmotic_pressure().magnitude - solution.get_osmotic_pressure().magnitude return (act_cation_mem/act_cation_soln) ** (1/z_cation) * (act_anion_soln/act_anion_mem)**(1/z_anion) - math.exp(delta_pi * exp_term) # solve the function above using one of scipy's nonlinear solvers from scipy.optimize import brentq # determine which ion concentration represents the co-ion # call a nonlinear solver to adjust the concentrations per the donnan # equilibrium, unless the membrane is uncharged # the initial guess is to set the co-ion concentration in the membrane # equal to that in the solution if fixed_charge.magnitude >0: x = conc_cation_soln.magnitude brentq(donnan_solve,1e-10,x,xtol=0.001) elif fixed_charge.magnitude <0: x = conc_anion_soln.magnitude brentq(donnan_solve,1e-10,x,xtol=0.001) else: pass # return the equilibrated solution return donnan_soln
def adjust_temp_vanthoff( equilibrium_constant, enthalpy, temperature, reference_temperature=25 * unit("degC") ): """(float,float,number, optional number) -> float Adjust a reaction equilibrium constant from one temperature to another. Parameters ---------- equilibrium_constant : float The reaction equilibrium constant for the reaction enthalpy : Quantity The enthalpy change (delta H) for the reaction in kJ/mol. Assumed independent of temperature (see Notes). temperature : Quantity the desired reaction temperature in degrees Celsius reference_temperature : Quantity, optional the temperature at which equilibrium_constant is valid. (25 degrees C if omitted). Returns ------- float adjusted reaction equilibrium constant Notes ----- This function implements the Van't Hoff equation to adjust measured equilibrium constants to other temperatures. .. math:: ln(K2 / K1) = {\\delta H \\over R} ( {1 \\over T_1} - {1 \\over T_2} ) This implementation assumes that the enthalpy is independent of temperature over the range of interest.[1] .. [1] Stumm, Werner and Morgan, James J. Aquatic Chemistry, 3rd ed, pp 53. Wiley Interscience, 1996. Examples -------- >>> adjust_temp_vanthoff(0.15,-197.6*unit('kJ/mol'),42*unit('degC'),25*unit('degC')) #doctest: +ELLIPSIS 0.00203566... If the 'ref_temperature' parameter is omitted, a default of 25 C is used. >>> adjust_temp_vanthoff(0.15,-197.6*unit('kJ/mol'),42*unit('degC')) #doctest: +ELLIPSIS 0.00203566... """ output = equilibrium_constant * math.exp( enthalpy / unit.R * (1 / reference_temperature.to("K") - 1 / temperature.to("K")) ) logger.info( "Adjusted equilibrium constant K=%s from %s to %s degrees Celsius with Delta H = %s. " "Adjusted K = %s % equilibrium_constant,reference_temperature,temperature,enthalpy,output" ) logger.warning( "Van't Hoff equation assumes enthalpy is independent of temperature over the range of interest" ) return output
def mix(Solution1, Solution2): ''' Mix two solutions together Returns a new Solution object that results from the mixing of Solution1 and Solution2 Parameters ---------- Solution1, Solution2 : Solution objects The two solutions to be mixed. Returns ------- Solution A Solution object representing the mixed solution. ''' # check to see if the two solutions have the same solvent if not Solution1.solvent_name == Solution2.solvent_name: logger.error('mix() function does not support solutions with different solvents. Aborting.') if not Solution1.solvent_name == 'H2O' or Solution1.solvent_name == 'water': logger.error('mix() function does not support non-water solvents. Aborting.') # set the pressure for the new solution p1 = Solution1.get_pressure() t1 = Solution1.get_temperature() v1 = Solution1.get_volume() p2 = Solution2.get_pressure() t2 = Solution2.get_temperature() v2 = Solution2.get_volume() # check to see if the solutions have the same temperature and pressure if not p1 == p2: logger.info('mix() function called between two solutions of different pressure. Pressures will be averaged (weighted by volume)') blend_pressure = str((p1 * v1 + p2 * v2) / (v1 + v2)) if not t1 == t2: logger.info('mix() function called between two solutions of different temperature. Temperatures will be averaged (weighted by volume)') blend_temperature = str((t1 * v1 + t2 * v2) / (v1 + v2)) # retrieve the amount of each component in the parent solution and # store in a list. mix_species={} for item in Solution1.components: mix_species.update({item:str(Solution1.get_amount(item,'mol'))}) for item in Solution2.components: if item in mix_species: new_amt = str(unit(mix_species[item]) + Solution2.get_amount(item,'mol')) mix_species.update({item:new_amt}) else: mix_species.update({item:Solution2.get_amount(item,'mol')}) # create an empty solution for the mixture Blend = pyEQL.Solution(temperature = blend_temperature,pressure= blend_pressure) # set or add the appropriate amount of all the components for item in mix_species.keys(): if item in Blend.components: # if already present (e.g. H2O, H+), modify the amount Blend.set_amount(item,mix_species[item]) else: # if not already present, add the component Blend.add_solute(item,mix_species[item]) return Blend
def adjust_temp_arrhenius( rate_constant, activation_energy, temperature, reference_temperature=25 * unit("degC"), ): """(float,float,number, optional number) -> float Adjust a reaction equilibrium constant from one temperature to another. Parameters ---------- rate_constant : Quantity The parameter value (usually a rate constant) being adjusted activation_energy : Quantity The activation energy of the process, in kJ/mol temperature : Quantity the desired reaction temperature. reference_temperature : Quantity, optional the temperature at which equilibrium_constant is valid Defaults to 25 degrees C if omitted. Returns ------- Quantity The adjusted reaction equilibrium constant See Also -------- kelvin Notes ----- This function implements the Arrhenius equation to adjust measured rate constants to other temperatures. [1] .. math:: ln(K2 / K1) = {E_a \\over R} ( {1 \\over T_1} - {1 \\over T_2} ) .. [1] http://chemwiki.ucdavis.edu/Physical_Chemistry/Kinetics/Reaction_Rates/\ Temperature_Dependence_of_Reaction_Rates/Arrhenius_Equation TODO - add better reference Examples -------- >>> adjust_temp_arrhenius(7,900*unit('kJ/mol'),37*unit('degC'),97*unit('degC')) #doctest: +ELLIPSIS 1.8867225...e-24 """ output = rate_constant * math.exp( activation_energy / unit.R * (1 / reference_temperature.to("K") - 1 / temperature.to("K")) ) logger.info( "Adjusted parameter %s from %s to %s degrees Celsius with Activation Energy = %s kJ/mol. " "Adjusted value = %s % rate_constant,reference_temperature,temperature,activation_energy,output" ) return output