pmutt.statmech.rot.RigidRotor
- class pmutt.statmech.rot.RigidRotor(symmetrynumber, rot_temperatures=None, geometry=None, atoms=None, degree_tol=5.0)
Bases:
_ModelBase
Rotational mode using the rigid rotor assumption. Equations sourced from:
Sandler, S. I. An Introduction to Applied Statistical Thermodynamics; John Wiley & Sons, 2010.
- symmetrynumber
Symmetry number of molecule. If a string is inputted, it represents the point group. Supported options are shown below.
Point group
symmetry number
C1
1
Cs
1
C2
2
C2v
2
C3v
3
Cinfv
1
D2h
4
D3h
6
D5h
10
Dinfh
2
D3d
6
Td
12
Oh
24
See DOI for more details: 10.1007/s00214-007-0328-0
- rot_temperatures
Rotational temperatures in K
- Type:
list of float, optional
- atoms
An atoms object can be used to calculate rot_temperatures and guess geometry
- Type:
ase.Atoms object, optional
- degree_tol
Degree tolerance to estimate geometry. Only required if estimating geometry or rot_temperatures
- Type:
- __init__(symmetrynumber, rot_temperatures=None, geometry=None, atoms=None, degree_tol=5.0)
Methods
__init__
(symmetrynumber[, rot_temperatures, ...])from_dict
(json_obj)Recreate an object from the JSON representation.
get_Cp
(units, **kwargs)Calculate the heat capacity (constant P)
get_CpoR
()Calculates the dimensionless heat capacity at constant pressure
get_Cv
(units, **kwargs)Calculate the heat capacity (constant V)
get_CvoR
()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
()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
()Calculates the dimensionless internal energy
get_q
(T)Calculates the partition function
to_dict
()Represents object as dictionary with JSON-accepted datatypes
- classmethod from_dict(json_obj)
Recreate an object from the JSON representation.
- Parameters:
json_obj (dict) – JSON representation
- Returns:
Obj
- Return type:
Appropriate object
- get_Cp(units, **kwargs)
Calculate the heat capacity (constant P)
- get_CpoR()
Calculates the dimensionless heat capacity at constant pressure
\(\frac{C_P^{rot}}{R}=\frac{C_V^{rot}}{R}\)
- Returns:
CpoR_rot – Rotational dimensionless heat capacity at constant pressure
- Return type:
- get_Cv(units, **kwargs)
Calculate the heat capacity (constant V)
- get_CvoR()
Calculates the dimensionless heat capacity at constant volume
\(\frac{C_V^{rot}}{R}=0\) if monatomic
\(\frac{C_V^{rot}}{R}=1\) if linear
\(\frac{C_V^{rot}}{R}=1.5\) if nonlinear
- Returns:
CvoR_rot – Rotational dimensionless heat capacity at constant volume
- Return type:
- get_F(units, T=298.15, **kwargs)
Calculate the Helmholtz energy
- get_FoRT(T)
Calculates the dimensionless Helmholtz energy
\(\frac{A^{rot}}{RT}=\frac{U^{rot}}{RT}-\frac{S^{rot}}{R}\)
- get_G(units, T=298.15, **kwargs)
Calculate the Gibbs energy
- get_GoRT(T)
Calculates the dimensionless Gibbs energy
\(\frac{G^{rot}}{RT}=\frac{H^{rot}}{RT}-\frac{S^{rot}}{R}\)
- get_H(units, T=298.15, **kwargs)
Calculate the enthalpy
- get_HoRT()
Calculates the dimensionless enthalpy
\(\frac{U^{rot}}{RT}=\frac{H^{rot}}{RT}\)
- Returns:
HoRT_rot – Rotational dimensionless enthalpy
- Return type:
- get_S(units, **kwargs)
Calculate the entropy
- get_SoR(T)
Calculates the dimensionless entropy
\(\frac{S^{rot}}{R}=0\) if monatomic
\(\frac{S}{R}=\log(\frac{T}{\sigma}\prod_i \frac{1}{ \Theta_{R,i}})+1\) if linear
\(\frac{S^{rot}}{R}=\log\bigg(\frac{\sqrt\pi }{\sigma} \prod_i{(\frac{T^3}{\Theta_{R,i}})^\frac{1}{2}}\bigg)+1.5\) if nonlinear
- get_U(units, T=298.15, **kwargs)
Calculate the internal energy
- get_UoRT()
Calculates the dimensionless internal energy
\(\frac{U^{rot}}{RT}=0\) if monatomic
\(\frac{U^{rot}}{RT}=1\) if linear
\(\frac{U^{rot}}{RT}=1.5\) if nonlinear
- Returns:
UoRT_rot – Rotational dimensionless internal energy
- Return type:
- get_q(T)
Calculates the partition function
\(q^{rot}=0\) if monatomic
\(q^{rot}=\frac{T}{\sigma}\prod_i\frac{1}{\Theta_R}\) if linear
\(q^{rot}=\frac{\sqrt{\pi}}{\sigma}(T^3\prod_i\frac{1}{ \Theta_{Ri}})^\frac{1}{2}\) if nonlinear