pmutt.statmech.vib.HarmonicVib
- class pmutt.statmech.vib.HarmonicVib(vib_wavenumbers=[], imaginary_substitute=None)
Bases:
_ModelBase
Vibrational modes using the harmonic approximation. Equations used sourced from:
Sandler, S. I. An Introduction to Applied Statistical Thermodynamics; John Wiley & Sons, 2010.
- vib_wavenumbers
Vibrational wavenumbers (\(\tilde{\nu}\)) in 1/cm
- Type:
list of float
- imaginary_substitute
If this value is set, imaginary frequencies are substituted with this value for calculations. Otherwise, imaginary frequencies are ignored. Default is None
- Type:
float, optional
- __init__(vib_wavenumbers=[], imaginary_substitute=None)
Methods
__init__
([vib_wavenumbers, imaginary_substitute])from_dict
(json_obj)Recreate an object from the JSON representation.
get_Cp
(units, **kwargs)Calculate the heat capacity (constant P)
get_CpoR
(T)Calculates the dimensionless heat capacity at constant pressure
get_Cv
(units, **kwargs)Calculate the heat capacity (constant V)
get_CvoR
(T)Calculates the dimensionless heat capacity at constant volume
get_F
(units[, T])Calculate the Helmholtz energy
get_FoRT
(T)Calculates the dimensionless Helmholtz energy
get_G
(units[, T])Calculate the Gibbs energy
get_GoRT
(T)Calculates the dimensionless Gibbs energy
get_H
(units[, T])Calculate the enthalpy
get_HoRT
(T)Calculates the dimensionless enthalpy
get_S
(units, **kwargs)Calculate the entropy
get_SoR
(T)Calculates the dimensionless entropy
get_U
(units[, T])Calculate the internal energy
get_UoRT
(T)Calculates the dimensionless internal energy
get_ZPE
()Calculates the zero point energy
get_q
(T[, include_ZPE])Calculates the partition function
Prints the wavenumbers that will be used in a thermodynamic calculation.
to_dict
()Represents object as dictionary with JSON-accepted datatypes
Attributes
- classmethod from_dict(json_obj)
Recreate an object from the JSON representation.
- Parameters:
json_obj (dict) – JSON representation
- Returns:
HarmonicVib
- Return type:
HarmonicVib object
- get_Cp(units, **kwargs)
Calculate the heat capacity (constant P)
- get_CpoR(T)
Calculates the dimensionless heat capacity at constant pressure
\(\frac{C_P^{vib}}{R}=\frac{C_V^{vib}}{R}=\sum_i \bigg(\frac{ \Theta_{V,i}}{2T}\bigg)^2 \frac{1}{\big(\sinh{\frac{\Theta_{V,i}} {2T}}\big)^2}\)
- get_Cv(units, **kwargs)
Calculate the heat capacity (constant V)
- get_CvoR(T)
Calculates the dimensionless heat capacity at constant volume
\(\frac{C_V^{vib}}{R}=\sum_i \bigg(\frac{\Theta_{V,i}}{2T} \bigg)^2 \frac{1}{\big(\sinh{\frac{\Theta_{V,i}}{2T}}\big)^2}\)
- get_F(units, T=298.15, **kwargs)
Calculate the Helmholtz energy
- get_FoRT(T)
Calculates the dimensionless Helmholtz energy
\(\frac{A^{vib}}{RT}=\frac{U^{vib}}{RT}-\frac{S^{vib}}{R}\)
- get_G(units, T=298.15, **kwargs)
Calculate the Gibbs energy
- get_GoRT(T)
Calculates the dimensionless Gibbs energy
\(\frac{G^{vib}}{RT}=\frac{H^{vib}}{RT}-\frac{S^{vib}}{R}\)
- get_H(units, T=298.15, **kwargs)
Calculate the enthalpy
- get_HoRT(T)
Calculates the dimensionless enthalpy
\(\frac{H^{vib}}{RT}=\frac{U^{vib}}{RT}=\sum_i \bigg(\frac{ \Theta_{V,i}}{2T}+\frac{\Theta_{V,i}}{T}\frac{\exp\big(-\frac{ \Theta_{V,i}}{T}\big)}{1-\exp\big(-\frac{\Theta_{V_i}}{T}\big)} \bigg)\)
- get_S(units, **kwargs)
Calculate the entropy
- get_SoR(T)
Calculates the dimensionless entropy
\(\frac{S^{vib}}{R}=\sum_i \frac{\Theta_{V,i}}{T}\frac{\exp \big(-\frac{\Theta_{V,i}}{T}\big)}{1-\exp\big(-\frac{ \Theta_{V,i}}{T}\big)}-\ln \bigg(1-\exp\big(-\frac{ \Theta_{V,i}}{T}\big)\bigg)\)
- get_U(units, T=298.15, **kwargs)
Calculate the internal energy
- get_UoRT(T)
Calculates the dimensionless internal energy
\(\frac{U^{vib}}{RT}=\sum_i \bigg(\frac{\Theta_{V,i}}{2T}+ \frac{\Theta_{V,i}}{T}\frac{\exp\big(-\frac{\Theta_{V,i}}{T} \big)}{1-\exp\big(-\frac{\Theta_{V_i}}{T}\big)}\bigg)\)
- get_ZPE()
Calculates the zero point energy
\(ZPE=\frac{1}{2}k_b\sum_i \Theta_{V,i}\)
- Returns:
zpe – Zero point energy in eV
- Return type:
- get_q(T, include_ZPE=True)
Calculates the partition function
\(q^{vib}=\prod_i \frac{\exp({-\frac{\Theta_{V,i}}{2T}})} {1-\exp({-\frac{\Theta_{V,i}}{T}})}\) if include_ZPE = True \(q^{vib}=\prod_i \frac{1} {1-\exp({-\frac{\Theta_{V,i}}{T}})}\) if include_ZPE = False
- print_calc_wavenumbers()
Prints the wavenumbers that will be used in a thermodynamic calculation. If
self.imaginary_substitute
is a float, then imaginary frequencies are replaced with that value. Otherwise, imaginary frequencies are ignored.
- to_dict()
Represents object as dictionary with JSON-accepted datatypes
- Returns:
obj_dict
- Return type:
- property vib_wavenumbers