Table Of Contents

Previous topic

pyviblib.calc.spectra

Next topic

pyviblib.gui

This Page

pyviblib.calc.vibrations

Module for handling vibrational motion

The formulae used to implement this modules have been taken from W. Hug and M. Fedorovsky, Theor. Chem. Acc. 119 (2006) 113-131, [ doi:10.1007/s00214-006-0185-2 ]

Functions

correlate_vibrations()
correlate two vibrations
remove_contaminations()
remove non-vibrational contaminations
correlate_vibrations_all()
correlate all found vibrations
create_dyad()
create the dyad for a vibration
generate_tr_rot()
generate translations/rotations
vibana()
perform vibrational analysis
omega_p()
angular frequency
dcm_p()
distance change matrix (DCM) for a vibration
vib_amplitudes()
amplitudes of atoms for a vibration
calc_all_devs()
conventional deviations of frequencies
hessian_dev()
deviation of the magnitudes of Hessians
hessian_magnitude()
magnitude of a Hessian
degree_localization()
degree of localization on a fragment
Author:Maxim Fedorovsky
pyviblib.calc.vibrations.correlate_vibrations(L_ref_p, L_tr_p, atom_pairs, U, remove_tr_rot=False, L_tr_rot_F_ref=None, L_tr_rot_F_tr=None)

Calculate the similarity and overlap of two vibrations on a fragment.

The vibrations must be given as mass-weighted excursions. Let us denote the first vibration as reference and the second one as trial.

Parameters:
  • L_ref_p – reference vibration (one-based ndarray), shape : (1 + Natoms_ref, 4) with Natoms_ref being the number of atoms in the reference molecule
  • L_tr_p – trial vibration (one-based ndarray), shape : (1 + Natoms_tr, 4) with Natoms_ref being the number of atoms in the trial molecule
  • atom_pairs – atom pairs comprising the fragment (null-based ndarray), shape : (Natoms_F, 2) with Natoms_F being the number of atoms in the fragment atom numbers are one-based, example : array([ [1, 1], [3, 3], [9, 9] ])
  • U – rotation matrix found by superimposing the reference and trial molecules on the fragment (one-based ndarray), shape : (4, 4), The trial vibration is rotated with this matrix. if None, the rotation is not carried out. refer to pyviblib.calc.qtrfit.fit()
  • remove_tr_rot – whether the non-vibrational contaminations are to be removed
  • L_tr_rot_F_ref – translations/rotations on the reference FRAGMENT (one-based ndarray), shape : (4 + nrot_ref, Natoms_F, 4) with nrot_ref being the number of rotational degrees of freedom in the reference molecule must be supplied if remove_tr_rot is True
  • L_tr_rot_F_tr – translations/rotations on the trial FRAGMENT (one-based ndarray) shape : (4 + nrot_tr, Natoms_F, 4) with nrot_tr being the number of rotational degrees of freedom in the trial molecule, must be supplied if remove_tr_rot is True
Returns:

dictionary with the following keys:

  • similarity - similarity
  • overlapoverlap
pyviblib.calc.vibrations.remove_contaminations(L_p, L_tr_rot_F, atoms)

Remove the translational and rotational contaminations on a fragment.

Parameters:
  • L_p – mass-weighted excursion for the vibration (one-based ndarray), shape : (1 + Natoms, 4) with Natoms being the number of atoms in the molecule
  • L_tr_rot_F – mass-weighted excursions for the translations/rotations on the fragment (one-based ndarray) shape : (4 + nrot, Natoms_F, 4) with nrot being the number of rotational degrees of freedom and Natoms_F the number of atoms in the fragment
  • atoms – atoms comprising the fragment (null-based ndarray), shape : (Natoms_F,), atom numbers are one-based
Returns:

the dyad for the vibration on the fragment, from which the rotational/translational contaminations have been removed, shape : (1 + Natoms_F, 4, 1 + Natoms_F, 4)

pyviblib.calc.vibrations.correlate_vibrations_all(L_ref, L_tr_rot_ref, L_tr, L_tr_rot_tr, atom_pairs, U=None, include_tr_rot=False, remove_tr_rot=False, L_tr_rot_F_ref=None, L_tr_rot_F_tr=None)

Correlate a set of vibrations on a fragment.

The vibrations to be correlated must be given as mass-weighted excursions. Let us denote by reference the first set of vibrations and by trial the second one.

Parameters:
  • L_ref – vibrations of the reference molecule (one-based ndarray), shape : (1 + NFreq_ref, 1 + Natoms_ref, 4) with NFreq_ref being the number of vibrations in the reference molecule and Natoms_ref the number of atoms in it
  • L_tr_rot_ref – translations/rotations of the reference molecule (one-based ndarray), shape : (4 + nrot_ref, 1 + Natoms_ref, 4) with nrot_ref being the number of rotational degrees of freedom in the reference molecule
  • L_tr – vibrations of the trial molecule (one-based ndarray), shape : (1 + NFreq_tr, 1 + Natoms_tr, 4) with NFreq_tr being the number of vibrations in the trial molecule and Natoms_tr the number of atoms in it
  • L_tr_rot_tr – translations/rotations of the trial molecule (one-based ndarray), shape : (4 + nrot_tr, 1 + Natoms_tr, 4) with nrot_tr being the number of rotational degrees of freedom in the trial molecule
  • atom_pairs – atom pairs comprising the fragment (null-based ndarray) shape : (Natoms_F, 2) with Natoms_F being the number of atoms in the fragment atom numbers are one-based example : array([ [1, 1], [3, 3], [9, 9] ])
  • U – rotation matrix found by superimposing the reference and trial fragments (one-based ndarray, default None), shape : (4, 4) The trial vibrations are rotated with this matrix. if None, the rotation is not carried out.
  • include_tr_rot – include the translations/rotations explicitely in the result overlaps and similarities matrices
  • remove_tr_rot – whether the non-vibrational contaminations are to be removed
  • L_tr_rot_F_ref – translations/rotations on the reference fragment (one-based ndarray) shape : (4 + nrot_ref, Natoms_F, 4) with nrot_ref being the number of rotational degrees of freedom in the reference molecule and Natoms_F the number of atoms in the fragment must be supplied if remove_tr_rot is True
  • L_tr_rot_F_tr – translations/rotations on the trial fragment (one-based ndarray), shape : (4 + nrot_tr, Natoms_F, 4) with nrot_tr being the number of rotational degrees of freedom in the trial molecule must be supplied if remove_tr_rot is True
Returns:

dictionary with the following keys:

  • overlaps - overlaps matrix (one-based ndarray)
  • similarities - similarities matrix (one-based ndarray)
pyviblib.calc.vibrations.create_dyad(L_p, atom_list_rows=None, atom_list_cols=None)

Create the dyad tensor for a vibration

Parameters:
  • L_p – mass-weighted excursion for the vibration (one-based ndarray), shape : (1 + Natoms, 4) with Natoms being the number of atoms in the molecule
  • atom_list_rows – atoms over which it is to be iterated through the rows (null-based ndarray), atom numbers are one-based, if None, use all atoms
  • atom_list_cols – atoms over which it is to be iterated through the columns (null-based ndarray), atom numbers are one-based, if None, use all atoms
Returns:

the dyad with the shape (1 + Natoms_rows, 4, 1 + Natoms_cols, 4)

pyviblib.calc.vibrations.generate_tr_rot(L, coords, masses, atom_list=None)

Generate translations/rotations for a fragment in a molecule.

Parameters:
  • L – set of vibrations (one-based ndarray), if given, the result translations/rotations will be made orthogonal to the set, shape : (1 + Nvib, 1 + Natoms, 4) with Nvib being number of vibrations, Natoms the number of atoms
  • coords – coordinates of the atoms (one-based ndarray), shape : (1 + Natoms, 4)
  • masses – masses of the atoms (one-based ndarray), shape : (1 + Natoms,)
  • atom_list – atoms comprising the fragment (null-based ndarray), atom numbers are one-based, shape : (Natoms_F,) with Natoms_F being the number of atoms in the fragment, if None, use all atoms

The return value is a dictionary with the following keys:

  • L_tr_rot - set of the translations and rotations (one-based ndarray)

    shape : (4 + nrot, 1 + Natoms_F, 4) with nrot being the number of rotational degrees of freedom. The first three members are always translations. The number of rotations nrot can vary.

  • I - principal axes of the inertia tensor (null-based ndarray),

    shape : (3,)

  • I_values - principal values of the inertia tensor

    (null-based ndarray), shape : (3,)

pyviblib.calc.vibrations.vibana(hessian, coords, masses)

Perform the vibrational analysis.

For details refer to http://www.gaussian.com/g_whitepap/vib.htm. This implementation generates the translations and rotations around the pricipal axes of the inertia tensor.

Parameters:
  • hessian – Hessian matrix (one-based ndarray), shape : (1 + Natoms, 4, 1 + Natoms, 4)
  • coords – coordinates of the atoms (one-based ndarray), shape : (1 + Natoms, 4)
  • masses – masses of the atoms (one-based ndarray), shape : (1 + Natoms,)
Returns:

dictionary with the following keys:

  • freqs - frequencies array sorted in ascending order

    (one-based ndarray), shape : (1 + NFreq,) with NFreq being the number of vibrations

  • L - set of mass-weighted excursions (one-based ndarray),

    shape : (1 + NFreq, 1 + Natoms, 4)

pyviblib.calc.vibrations.omega_p(nu)

Calculate the angular frequency.

\omega_p = 200 \pi c \tilde{\nu}

Parameter:nu – wavenumber in cm**(-1)

The return value is expressed in inverse seconds.

pyviblib.calc.vibrations.dcm_p(coords, Lx, nu)

Calculate the distance change matrix (DCM) for a normal mode.

Parameters:
  • coords – coordinates of the atoms (one-based ndarray), shape : (1 + Natoms, 4) with Natoms being the number of atoms
  • Lx – cartesian excursions (one-based ndarray), shape : (1 + Natoms, 4)
  • nu – wavenumber in cm**(-1)
Returns:

one-based ndarray with the shape (1 + Natoms, 1 + Natoms)

pyviblib.calc.vibrations.vib_amplitudes(Lx, nu, n_p=0)

Calculate the amplitudes of atoms for a vibration.

Parameters:
  • Lx – cartesian excursions (one-based ndarray), shape : (1 + Natoms, 4) with Natoms being the number of atoms
  • nu – wavenumber in cm**(-1)
  • n_p – vibrational number, the default value instructs to calculate the zero-point amplitudes
Returns:

one-based ndarray with the shape (1 + Natoms,), the units are angstroms

pyviblib.calc.vibrations.calc_all_devs(ref_freqs, tr_freqs)

Calculate deviations of trial frequencies from reference ones.

Parameters:
  • ref_freqs – reference frequencies (one-dimensional ndarray)
  • tr_freqs – trial frequencies (one-dimensional ndarray)
The following magnitudes are calculated:
  • ALG2 - relative deviation of the algebraic mean of the quadrats of

    the frequencies (%)

  • GEOM2 - relative deviation of the geometric mean of the quadrats of

    the frequencies (%)

  • MD - mean deviation of the frequencies (units of the frequencies)

  • MAD - mean absolute deviation of the frequencies

    (units of the frequencies)

  • MAXD - maximal absolute deviation of the frequencies

    (units of the frequencies)

  • RELMD - mean relative deviation of the frequencies (%)

  • RELMAD - mean relative absolute deviation of the frequencies (%)

  • RELMAXD - maximal relative absolute deviation of the frequencies (%)

  • FRMSD - root mean square deviations of the frequencies (%)

Returns:(ALG2, GEOM2, MD, MAD, MAXD, RELMD, RELMAD, RELMAXD, FRMSD)
pyviblib.calc.vibrations.hessian_dev(ref_f, tr_f, row_atoms=None, col_atoms=None, rotmat=None)

Calculate the normalized deviation of the magnitudes of Hessians.

If a rotation matrix is supplied then the trial Hessian is to be rotated with it.

Parameters:
  • ref_f – reference Hessian (one-based ndarray)
  • tr_f – trial Hessian (one-based ndarray)
  • row_atoms – one-based list of row atoms (iteratable array) if None, all atoms are used
  • col_atoms – one-based list of column atoms (iteratable array) if None, all atoms are used
  • rotmat – left rotation matrix (one-based ndarray), shape : (4, 4)
pyviblib.calc.vibrations.hessian_magnitude(f, row_atoms=None, col_atoms=None)

Calculate the magnitude of a Hessian.

Parameters:
  • f – Hessian (one-based ndarray)
  • row_atoms – one-based list of row atoms (iteratable array), if None, all atoms are used
  • col_atoms – one-based list of column atoms (iteratable array), if None, all atoms are used
pyviblib.calc.vibrations.degree_localization(L_p, atom_list)

Degree of localization of a vibration on a fragment.

Parameters:
  • L_p – mass-weighted excursion for the vibration (one-based ndarray), shape : (1 + Natoms, 4) with Natoms being the number of atoms in the molecule
  • atom_list – atoms comprising the fragment (null-based ndarray), atom numbers are one-based shape : (Natoms_F,) with Natoms_F being the number of atoms in the fragment
Returns:

value within [0, 1]


SourceForge.net Logo