vibrational_analysis

Vibrational analysis.

class daltonproject.vibrational_analysis.VibrationalAnalysis(frequencies: np.ndarray, ir_intensities: np.ndarray | None = None, raman_intensities: np.ndarray | None = None, cartesian_displacements: np.ndarray | None = None)

Data structure for vibrational properties.

cartesian_displacements: np.ndarray | None

Alias for field number 3

frequencies: np.ndarray

Alias for field number 0

ir_intensities: np.ndarray | None

Alias for field number 1

raman_intensities: np.ndarray | None

Alias for field number 2

daltonproject.vibrational_analysis.compute_ir_intensities(dipole_gradients: ndarray, transformation_matrix: ndarray, is_linear: bool = False) ndarray

Compute infrared molar decadic absorption coefficient in m^2 / (s * mol).

Parameters
  • dipole_gradients – Gradients of the dipole moment (au) with respect to nuclear displacements in Cartesian coordinates.

  • transformation_matrix – Transformation matrix to convert from Cartesian to normal coordinates.

  • is_linear – Indicate if the molecule is linear.

Returns

IR intensities in m^2 / (s * mol).

daltonproject.vibrational_analysis.compute_raman_intensities(polarizability_gradients: ndarray, polarizability_frequency: float, vibrational_frequencies: ndarray, transformation_matrix: ndarray, is_linear: bool = False) ndarray

Compute Raman differential scattering cross-section in C^4 * s^2 / (J * m^2 * kg).

Parameters
  • polarizability_gradients – Gradients of the polarizability with respect to nuclear displacements in Cartesian coordinates.

  • polarizability_frequency – Frequency of the incident light in au.

  • vibrational_frequencies – Vibrational frequencies in cm^-1.

  • transformation_matrix – Transformation matrix to convert from Cartesian to normal coordinates.

  • is_linear – Indicate if the molecule is linear.

Returns

Raman intensities in C^4 * s^2 / (J * m^2 * kg).

daltonproject.vibrational_analysis.normal_coords(molecule: Molecule, hessian: np.ndarray, is_linear: bool = False, isotopes: Sequence[float] | None = None) np.ndarray

Compute normal coordinates.

Parameters
  • molecule – Molecule object.

  • hessian – Second derivatives of energy with respect to nuclear displacements.

  • is_linear – Indicate if the molecule is linear.

  • isotopes – List of isotopes for each atom in the molecule.

Returns

q_mat = normal mode coordinates.

daltonproject.vibrational_analysis.normal_mode_eigensolver(molecule: Molecule, hessian: np.ndarray, is_linear: bool = False, isotopes: Sequence[float] | None = None) NormalModeEigenSolution

Calculate eigenvalues and eigenvectors of the Hessian matrix.

Parameters
  • molecule – Molecule object.

  • hessian – Second derivatives of energy with respect to nuclear displacements in Cartesian coordinates.

  • is_linear – Indicate if the molecule is linear.

  • isotopes – List of isotopes for each atom in the molecule.

Returns

Eigenvalues of the mass-weighted Hessian matrix. transformation_matrix: Eigenvectors of the mass-weighted Hessian matrix. frequencies: Vibrational frequencies in cm^-1.

Return type

eigen_vals

daltonproject.vibrational_analysis.number_modes(molecule: Molecule, is_linear: bool) int

Return number of normal modes.

Parameters
  • molecule – Molecule object.

  • is_linear – Indicate if the molecule is linear.

Returns

Number of normal modes.

Return type

num_modes

daltonproject.vibrational_analysis.vibrational_analysis(molecule: Molecule, hessian: np.ndarray, dipole_gradients: np.ndarray | None = None, polarizability_gradients: PolarizabilityGradients | None = None) VibrationalAnalysis

Perform vibrational analysis.

Parameters
  • molecule – Molecule object.

  • hessian – Second derivatives with respect to nuclear displacements in Cartesian coordinates.

  • dipole_gradients – Gradient of the dipole moment with respect to nuclear displacements in Cartesian coordinates.

  • polarizability_gradients – Gradient of the polarizability with respect to nuclear displacements in Cartesian coordinates.

Returns

Harmonic vibrational frequencies (cm^-1). Optionally, IR intensities (m^2 / (s * mol)), Raman intensities (C^4 * s^2 / (J * m^2 * kg)), and transformation matrix.