Table Of Contents

Previous topic

pyviblib.calc.qtrfit

Next topic

pyviblib.calc.vibrations

This Page

pyviblib.calc.spectra

Module for Raman/ROA as well as IR/VCD spectra generation

Inheritance diagram of pyviblib.calc.spectra

The formulae used to implement this module have been taken from [1].

References

[1]W. Hug. Chem. Phys., 264(1):53-69, 2001, [ doi: 10.1016/S0301-0104(00)00390-6 ]
[2](1, 2) W. Hug and J. Haesler, Int. J. Quant. Chem. 104:695-715, 2005, [ doi: 10.1002/qua.20600 ]
[3]Atkins, 6th edition, p. 459.

Classes

RamanROATensors
manipulating the Raman/ROA tensors
IRVCDTensors
manipulating the IR/VCD tensors
ExpRamanROAData
experimental Raman, ROA or Degree of circularity spectra

Functions

contract_A_tensor()
contract the derivatives of the A tensor
f_Qp_i()
transition integral between two states
Kp()
constant Kp used for the Raman/ROA cross-sections
conv_units_raman_roa()
convenience function for conversions
conv_invariant()
conversion: a.u. -> invariants’ units
conv_cross_sections()
conversion: a.u. -> cross-section’ units
conv_A4_AMU()
conversion: a.u. -> A^4 / AMU
inv_coeffs_raman_roa()
linear combinations for the cross-sections
stick_raman_roa()
generate Raman/ROA stick spectra
fit_raman_roa()
generate Raman/ROA spectra as curves
make_acp()
generate atomic contribution patterns (ACPs)
calc_spectra_nmol()
calculate VROA spectra for several molecules

Module variables

LIST_INVARIANTS
list of available invariants
LIST_VTENSORS
list of available V tensors
LIST_INVARIANTSUNITS
list of available units of the invariants (J)
X_PEAK_INTERVAL
length of a peak’s fitting interval
X_STEP
default x step for fitting curves
Author:Maxim Fedorovsky
class pyviblib.calc.spectra.RamanROATensors(PP, PM, PQ, Lx, freqs, lambda_incident, A=None, gauge_origin=None, London_orbitals=True, molecule_id=None, coords=None)

Class for manipulating the Raman/ROA tensors.

Refer to the DALTON manual page for information how to request a ROA job. http://www.kjemi.uio.no/software/dalton/resources/dalton20manual.pdf

The following read-only properties are exposed:
Natoms
number of atoms
NFreq
number of vibrations
Lx
cartesian excursions for the vibrations
freqs
wavenumbers of vibrations in ascending order
PP
gradients of the polarizability tensor
PM
gradients of the G’ tensor
PQ
gradients of the A tensor (contracted)
A
gradients of the A tensor (non-contracted)
lambda_incident
wavelength of the incident light
gauge_origin
gauge origin
London_orbitals
whether the London orbitals are used
molecule_id
associated molecule_id
The following public methods are exported:
V_a2()
V-tensor for the a^2 molecular invariant
V_b2()
V-tensor for the b^2 molecular invariant
V_aG()
V-tensor for the aG' molecular invariant
V_b2G()
V-tensor for the b^2_G molecular invariant
V_b2A()
V-tensor for the b^2_A molecular invariant
V_spacial_odc()
Spacial term for the b^2_G and b^2_A invariants
set_V_tensors()
overwrite the V-tensors
a2()
matrix of the a^2 molecular invariant
b2()
matrix of the b^2 molecular invariant
aG()
matrix of the aG' molecular invariant
b2G()
matrix of the b^2_G molecular invariant
b2A()
matrix of the b^2_A molecular invariant
J_all()
full five molecular invariants for all vibrations
J_all_AB()
ditto but for two groups
J_all_dec()
decomposed matrices of the invariants
J_rot()
molecular invariants for the rotations
acp()
generate atomic contribution patterns (ACPs)
gcp()
generate group contribution patterns (GCPs)
intensities()
calculate the Raman/ROA intensities
fit_spectra()
calculate fitted spectra
stick_spectra()
calculate stick spectra
norm_sum_vib()
normalized sum of ROA over the vibrations
norm_sum_rot()
normalized sum of ROA over the rotations
make_enantiomer()
make the virtual enantiomer
transform()
perform a geometrical transformation
The following public static method is exported:
V_all()
calculate the V-tensors from the gradients
J_all_mol()
invariants for the vibrations of a molecule

Initializer of the class.

PP, Lx, freqs and lambda_incident are required to be set. If PM or PQ are not set, they will be filled with zeros.

Parameters:
  • PP – gradients of the polarizability tensor (one-based ndarray), shape: (1 + Natoms, 4, 4, 4) with Natoms being the number of atoms
  • PM – gradients of the G’ tensor (one-based ndarray), shape: (1 + Natoms, 4, 4, 4)
  • PQ – gradients of the A tensor (contracted) (one-based ndarray), shape: (1 + Natoms, 4, 4, 4)
  • Lx – cartesian excursions for the vibrations (one-based ndarray) shape: (1 + NFreq, 1 + Natoms, 4), NFreq is the number of vibrations
  • freqs – wavenumbers of vibrations in ascending order (one-based ndarray), shape: (1 + NFreq,)
  • lambda_incident – wavelength of the incident light in nm
  • A – gradients of the non-contracted A tensor (one-based ndarray, default None), shape: (1 + Natoms, 4, 4, 4, 4)
  • gauge_origin – gauge origin (ndarray, shape: (4,))
  • London_orbitals – whether the London orbitals are used
  • molecule_id – associated molecule_id
  • coords – cartesian coordinates (one-based ndarray) shape: (1 + Natoms, 4)
J_all(*args, **kw)

Calculate the five molecular invariants for all vibrations.

A row is made up of a^2, b^2, aG', b^2_G and b^2_A.

Parameters:
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
  • parallel – whether to use parallel calculation
  • distance_terms – distance terms (for b2G and b2A) * 0 – usual (local term + non-local term) * 1 – local term * 2 – non-local term
Returns:

zero-based ndarray with a shape of (NFreq, 5) with NFreq being the number of vibrations.

J_all_AB(*args, **kw)

Calculate the five molecular invariants for all vibrations taking into account the contributions of groups A and B.

A row is made up of a^2, b^2, aG', b^2_G and b^2_A.

Parameters:
  • group_A – ndarray with one-based numbers of atoms of group A
  • group_B – ndarray with one-based numbers of atoms of group B
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
  • parallel – whether to use parallel calculation
  • distance_terms – distance terms (for b2G and b2A) * 0 – usual (local term + non-local term) * 1 – local term * 2 – non-local term
Returns:

zero-based ndarray with a shape of (3, NFreq, 5) with NFreq being the number of vibrations. The first index refers to the contribution of AA (0), AB(1) and BB(2).

J_all_dec(parts, units='invariant', parallel=False)

Decompose the five molecular invariants for all vibrations.

A row is made up of a^2, b^2, aG', b^2_G and b^2_A with a shape of (nparts, nparts) where nparts is the number of parts.

Parameters:
  • parts – tuple or list of ndarrays of one-based numbers of atoms, which comprise the parts
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
  • parallel – whether to use parallel calculation
Returns:

ndarray with a shape of (NFreq, 5, nparts, nparts) with NFreq being the number of vibrations.

static J_all_mol(mol, units='cross-section', parallel=False)

Get the five molecular invariants for a given molecule.

Parameters:
  • mol – molecule
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
  • parallel – whether to use parallel calculation
J_rot(Lx_rot, parallel=False)

Calculate the five reduced molecular invariants for the rotations.

A row is made up of a^2, b^2, aG', b^2_G and b^2_A.

Parameters:
  • Lx_rot – rotations, shape: (nrot, 1 + Natoms, 4)
  • parallel – whether to use parallel calculation

Lx_rot can be obtained as follows:

Lx_rot = mol.Lx_tr_rot[1 + mol.nrot:]
Returns:zero-based ndarray with the shape of (nrot, 5) with nrot being the number of rotations.
V_a2(i=0, distance_terms=0)

Calculate the V-tensor for the a^2 molecular invariant.

Parameters:
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • distance_terms – distance terms * 0 – usual * 1 – local term, the same as usual * 2 – non-local term, zero per definition

The tensor is expressed in a.u.

V_aG(i=0, distance_terms=0)

Calculate the V-tensor for the aG' molecular invariant.

Parameters:
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • distance_terms – distance terms * 0 – usual * 1 – local term, the same as usual * 2 – non-local term, zero per definition

The tensor is expressed in a.u.

static V_all(PP, PM, PQ, lambda_incident, row_atoms=None, column_atoms=None, parallel=False)

Calculate the five V-tensors from passed gradients.

Parameters:
  • row_atoms – list or tuple of one-based number of atoms through which it is to be iterated, unless given, use all atoms
  • column_atoms – dito for the columns
  • parallel – whether to use parallel calculation
Returns:

V_a2, V_b2, V_aG, V_b2G, V_b2A

V_b2(i=0, distance_terms=0)

Calculate the V-tensor for the b^2 molecular invariant.

Parameters:
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • distance_terms – distance terms * 0 – usual * 1 – local term, the same as usual * 2 – non-local term, zero per definition

The tensor is expressed in a.u.

V_b2A(i=0, distance_terms=0)

Calculate the V-tensor for the b^2_A molecular invariant.

Parameters:
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • distance_terms – distance terms (for b2G and b2A) * 0 – usual (local term + non-local term) * 1 – local term * 2 – non-local term

The tensor is expressed in a.u.

V_b2G(i=0, distance_terms=0)

Calculate the V-tensor for the b^2_G molecular invariant.

Parameters:
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • distance_terms – distance terms (for b2G and b2A) * 0 – usual (local term + non-local term) * 1 – local term * 2 – non-local term

The tensor is expressed in a.u.

V_spacial_odc(i=0)

Calculate the spacial term for b^2_G and b^2_A.

Parameter:i

part of the decomposed tensor (default 0, i.e. total sum) possible values:

  • 0: total sum
  • 1: isotropic part
  • 2: anisotropic part
  • 3: antisymmetric part

The tensor is expressed in a.u.

a2(p, i=0, units='invariant')

Calculate the a^2 molecular invariant.

Parameters:
  • p – number of vibration
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
aG(p, i=0, units='invariant')

Calculate the aG' molecular invariant.

Parameters:
  • p – number of vibration
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
acp(item, p)

Generate atomic contribution patterns (ACPs).

See also make_acp().

Parameters:
  • item – moleculular invariant or a cross-section, one of (‘raman’, ‘roa_backward’, ‘roa_forward’, ‘a2’, ‘b2’, ‘aG’, ‘b2G’, ‘b2A’)
  • p – number of vibration

The ACPs have a unit of cross-sections for raman, roa_* etc. and that of invariants for single invariants.

b2(p, i=0, units='invariant')

Calculate the b^2 molecular invariant.

Parameters:
  • p – number of vibration
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
b2A(p, i=0, units='invariant')

Calculate the b^2_A molecular invariant.

Parameters:
  • p – number of vibration
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
b2G(p, i=0, units='invariant')

Calculate the b^2_G molecular invariant.

Parameters:
  • p – number of vibration
  • i

    part of the decomposed tensor (default 0, i.e. total sum) possible values:

    • 0: total sum
    • 1: isotropic part
    • 2: anisotropic part
    • 3: antisymmetric part
  • units – units of the V-tensor, one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’)
fit_spectra(*args, **kw)

Calculate fitted spectra.

Parameters:
  • scattering – scattering type (see inv_coeffs_raman_roa())
  • N_G – number of Gauss functions
  • FWHM_is – full width at half maximum for the isotropic invariants a^2 and b^2
  • FWHM_anis – full width at half maximum for the anisotropic invariants aG', b^2_G, b^2_A
  • FWHM_inst – full width at half maximum for the Gaussian instrument profile
  • startx – start wavenumber of the fitting interval
  • endx – end wavenumbers of the fitting interval if the latter two keyword arguments are None, they are generated automatically
  • stepx – step in cm**(-1), the internal default value is used unless given
  • X – use these x values instead of generating it
  • ab – if not None, scale the frequencies as follows:
a = ab[0]
b = ab[1]
new_freqs = (a * freqs + b) * freqs
Parameters:
  • groups_AB – unless None, defines two groups A and B whose contributions are to be used for calculating the spectra, format : [group_A, group_B, t] t = 0 for AA, 1 for AB + BA and 2 for BB
  • contr_raman

    defines how the Raman intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of a2
    • 2 – only the contribution of b2
  • contr_roa

    defines how the ROA intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of aG’
    • 2 – only the contribution of b2G
    • 3 – only the contribution of b2A
  • contr_weight – whether to weight the Raman and ROA intensities. Only relevant if 0 < contr_raman or 0 < contr_roa
  • distance_terms – distance terms (for b2G and b2A) * 0 – usual (local term + non-local term) * 1 – local term * 2 – non-local term
Returns:

X, raman_Y, roa_Y, degcirc_Y

gcp(groups, item, p)

Generate group contribution patterns (GCPs).

Parameters:
  • groups – groups (list), atom indices are one-based, example : [ [1, 2], [4, 5], [6, 3] ]
  • item – moleculular invariant or a cross-section, one of (‘raman’, ‘roa_backward’, ‘roa_forward’, ‘a2’, ‘b2’, ‘aG’, ‘b2G’, ‘b2A’)
  • p – number of vibration

The GCPs have a unit of cross-sections for raman, roa_* etc. and that of invariants for single invariants.

static intensities(invars, scattering)

Calculate the Raman/ROA intensities.

Parameters:
  • invars – array with the five invariants (scalars or matrices)
  • scattering – scattering type (see inv_coeffs_raman_roa())
make_enantiomer()
Make the virtual enantiomer.
norm_sum_rot(masses, scattering=0)

Calculate the normalized sum of ROA over the rotations.

See eq 33 in [2].

Parameters:
  • masses – masses of the atoms in AMU (one-based ndarray), shape: (1 + Natoms,)
  • scattering – scattering type (see inv_coeffs_raman_roa())
norm_sum_vib(scattering=0, start_p=1, end_p=None)

Calculate the normalized sum of ROA over a range of vibrations.

It is a measure of the influence of the chiral distribution of a molecule’s electrons on vibrational ROA alone.

See eq 34 in [2].

Parameters:
  • scattering – scattering type (see inv_coeffs_raman_roa())
  • start_p – start vibration
  • end_p – end vibration, use None to specify the last
Returns:

freqs, roa_cumsum having a length of end_p - start_p + 1

set_V_tensors(V_tensors)

Overwrite the V-tensors generated from the gradients.

Parameter:V_tensors – list or tuple with V_a2, V_b2, V_aG, V_b2G, V_b2A, each of which is a ndarray with a shape of (1 + Natoms, 4, 1 + Natoms, 4)

An exception is raised unless the V-tensors are valid.

stick_spectra(scattering=0, groups_AB=None, contr_raman=0, contr_roa=0, contr_weight=True, distance_terms=0)

Calculate stick spectra.

Parameters:
  • scattering – scattering type (see inv_coeffs_raman_roa())
  • groups_AB – unless None, defines two groups A and B whose contributions are to be used for calculating the spectra, format : [group_A, group_B, t] t = 0 for AA, 1 for AB + BA and 2 for BB
  • contr_raman

    defines how the Raman intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of a2
    • 2 – only the contribution of b2
  • contr_roa

    defines how the ROA intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of aG’
    • 2 – only the contribution of b2G
    • 3 – only the contribution of b2A
  • contr_weight – whether to weight the Raman and ROA intensities. Only relevant if 0 < contr_raman or 0 < contr_roa
  • distance_terms – distance terms (for b2G and b2A) * 0 – usual (local term + non-local term) * 1 – local term * 2 – non-local term
Returns:

X, raman_Y, roa_Y, degcirc_Y. Each of these arrays is zero-based and of the length NFreq.

transform(rotmat, new_origin=None, trans_grad=False)

Perform a geometrical transformation.

New origin must be given if a translation of the origin dependent gradients of the G’ and A tensors is desired.

Parameters:
  • rotmat (ndarray, shape: (4, 4)) – rotation matrix
  • new_origin (ndarray, shape: (4,)) – new origin (in Angstrom)
  • tran_grad – whether to translate the gradients

The coordinates are assumed to be given in angstroms.

class pyviblib.calc.spectra.IRVCDTensors(P, M, Lx, freqs, molecule_id=None)

Class for manipulating the IR/VCD tensors.

Refer to the DALTON manual page for information on how to request a VCD job, http://www.kjemi.uio.no/software/dalton/resources/dalton20manual.pdf

The following read-only properties are exposed:
Natoms
number of atoms
NFreq
number of vibrations
Lx
cartesian excursions for the vibrations
freqs
wavenumbers of vibrations in ascending order
P
gradients of the electric dipole moment (APTs)
M
gradients of the magnetic dipole moment (AATs)
V_D
dyadics for the dipole strength
V_R
dyadics for the rotational strength
molecule_id
associated molecule_id
The following public methods are exported:
D()
matrix of the reduced dipole strength
R()
matrix of the reduced rotational strength
ir_matrix()
matrix of the infrared absorption intensity
vcd_matrix()
matrix of the VCD intensity
integrated_coeffs()
integrated absorbtion coefficients and g-factor
acp()
generate atomic contribution patterns (ACPs)
fit_spectra()
generate the IR/VCD/g in the curve repr.
norm_sum_vib()
normalized sum of VCD over the vibrations
norm_sum_rot()
normalized sum of VCD over the rotations
make_enantiomer()
make the virtual enantiomer
rotate()
perform a rotation

Initializer of the class.

P, Lx and freqs are required to be set. If M is not set, it will be filled with zeros.

Parameters:
  • P – gradients of the electric dipole moment (one-based ndarray), known in the VCD literature also as atomic polar tensors (APTs), shape: (1 + Natoms, 4, 4), Natoms is the number of atoms
  • M – gradients of the magnetic dipole moment (one-based ndarray) known in the VCD literature also as atomic axial tensors (AATs), shape: (1 + Natoms, 4, 4), Natoms is the number of atoms
  • Lx – cartesian excursions for the vibrations (one-based ndarray), shape: (1 + NFreq, 1 + Natoms, 4), NFreq is the number of vibrations
  • freqs – wavenumbers of vibrations in ascending order (one-based ndarray), shape: (1 + NFreq,)
  • molecule_id – associated molecule_id
D(p)

Calculate the matrix of the reduced dipole strength.

Parameter:p – number of vibration
R(p)

Calculate the matrix of the reduced rotational strength.

Parameter:p – number of vibration
acp(item, p)

Generate atomic contribution patterns (ACPs).

Parameters:
  • item – one of (‘IR’, ‘VCD’)
  • p – number of vibration

See also

make_acp()

fit_spectra(*args, **kw)

Generate the IR/VCD spectra as curves.

Parameters:
  • ngauss – number of Gauss functions
  • fwhm_fit – full width at half maximum of the Lorentz bands
  • fwhm_inst – full width at half maximum of the Gaussian instrument profile
  • startx – start wavenumber of the fitting interval
  • endx – end wavenumbers of the fitting interval, if the latter two keyword arguments are None, they are generated automatically
  • stepx – step in cm**(-1), the internal default value is used unless given
  • X – use these x values instead of generating it
  • ab – if not None, scale the frequencies as follows:
a = ab[0]
b = ab[1]
new_freqs = (a * freqs + b) * freqs
Parameter:groups_AB – unless None, defines two groups A and B whose contributions are to be used for calculating the spectra, format : [group_A, group_B, t] t = 0 for AA, 1 for AB + BA and 2 for BB
Returns:X, eps_Y, deps_Y, g_Y
integrated_coeffs(groups_AB=None)

Get the integrated absorption coeffcients (in SI) and the g-factors.

The shape of the bands is considered to be Lorentzian.

Parameter:groups_AB – unless None, defines two groups A and B whose contributions are to be used for calculating the spectra, format : [group_A, group_B, t] t = 0 for AA, 1 for AB + BA and 2 for BB

See [3].

Returns:ndarray with a shape of (NFreq, 3)
ir_matrix(p)

Calculate the matrix of the infrared absorption intensity.

Parameter:p – number of vibration

The units are m/mol.

make_enantiomer()
Make the virtual enantiomer.
norm_sum_rot(masses)

Calculate the normalized sum of VCD over the rotations.

Parameter:masses – masses of the atoms in AMU (one-based ndarray), shape: (1 + Natoms,)
norm_sum_vib(start_p=1, end_p=None)

Calculate the normalized sum of VCD over a range of vibrations.

Parameters:
  • start_p – start vibration
  • end_p – end vibration, use None to specify the last vibration
Returns:

freqs, roa_cumsum having a length of end_p - start_p + 1

rotate(rotmat)

Perform a rotation.

Parameter:rotmat (ndarray, shape: (4, 4)) – rotation matrix
vcd_matrix(p)

Calculate the matrix of the VCD intensity.

Parameter:p – number of vibration

The units are m/mol.

class pyviblib.calc.spectra.ExpRamanROAData(filename, laser_power, exposure_time, datatype, scattering)

Storing experimental Raman, ROA or Degree of circularity spectra.

The following read-only properties are exposed:
filename
file name
laser_power
laser power (mW)
exposure_time
exposure time (min)
datatype
see pyviblib.gui.esources.STRINGS_EXP_RAMANROA_SPECTRA_TYPES
scattering
see pyviblib.gui.resources.STRINGS_SCATTERING_TYPES
wavenumbers
wavenumbers
raman
raman
roa
roa
degcirc
degree of circularity

Initializer of the class.

Parameters:
  • filename – file name
  • laser_power – laser power (mW)
  • exposure_time – exposure time (min)
  • datatype – datatype
  • scattering – scattering type
pyviblib.calc.spectra.contract_A_tensor(A)

Contract the derivaritives of the A tensor.

The contraction with antisymmetric unit tensor of Levi-Civita \epsilon_{\mu\rho\sigma} is performed using the following formula (the Einstein summing convention of indices is used):

\frac{\partial \mathcal{A}^a_{\mu\nu}}{\partial x^j} =
\epsilon_{\mu\rho\sigma}
\frac{\partial A^a_{\rho\sigma\nu}}{\partial x^j}

Parameter:A – derivaritives of the A tensor (one-based ndarray), shape: (1 + Natoms, 4, 4, 4, 4) with Natoms being the number of atoms
Returns:the derivatives of the \mathcal{A} tensor with a shape of (1 + Natoms, 4, 4, 4)
pyviblib.calc.spectra.f_Qp_i(nu)

Transition integral between two vibrational states |i> and <f|.

\left\langle f | Q_p | i \right\rangle =
    \sqrt{\frac{\hbar}{2 \omega_p}}

It is assumed that the f <- i is a fundamental transition and that the vibration is harmonic.

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

Return value is expressed in SI units.

pyviblib.calc.spectra.Kp(nu_0, nu)

Factor appearing in the calculating of the Raman/ROA cross-sections.

K_p = \frac{10^7}{9} \pi^2 \mu_0^2 c^4 \left(\tilde{\nu}_0 -
    \tilde{\nu} \right)^3 \tilde{\nu}_0

Parameters:
  • nu_0 – wavenumber of the incident light in cm**(-1)
  • nu – shift in cm**(-1)

Return value is expressed in SI units.

pyviblib.calc.spectra.conv_units_raman_roa(units, lambda_incident, nu, J)

Conversion factor for a Raman/ROA quantity.

Parameters:
  • units – one of (‘invariant’, ‘cross-section’, ‘A^4/AMU’, ‘reduced’)
  • lambda_incident – wavelength of the incident light in nm
  • nu – wavenumber in cm**(-1)
  • J – molecular invariant (if applicable), one of (‘a2’, ‘b2’, ‘aG’, ‘b2G’, ‘b2A’)
pyviblib.calc.spectra.conv_invariant(nu)

Conversion factor from J_ab in a.u. to molecular invariants in SI units.

Parameter:nu – wavenumber in cm**(-1)
pyviblib.calc.spectra.conv_cross_sections(lambda_incident, nu)

Conversion factor from J_ab in a.u. to the cross-sections in A^2 / sr.

Parameters:
  • lambda_incident – wavelength of the incident light in nm
  • nu – wavenumber in cm**(-1)
pyviblib.calc.spectra.conv_A4_AMU(J)

Conversion factor from J_ab in a.u. to A^4/AMU.

For a2 and b2 (Raman invariants) : 142.943570. For aG, b2G, b2A (ROA invariants): 142.943570 / C_AU, with C_AU being the speed of light in a.u.

Parameter:J – molecular invariant, one of (‘a2’, ‘b2’, ‘aG’, ‘b2G’, ‘b2A’)
pyviblib.calc.spectra.inv_coeffs_raman_roa(scattering=0)

Linear combination coefficients of the five molecular invariants for the Raman/ROA cross-sections for a given scattering type.

\sigma & = K_{raman} \left(x_{1ram} a^2 + x_{2ram} b^2 \right) \\
\Delta \sigma & = K_{roa} \left(x_{1roa} aG + x_{2roa} b^2_G +
    x_{3roa} b^2_A\right)

Parameter:scattering – scattering, the default one is backward

Possible values of the scattering argument:

  • 0: backward ( 1., 90., 14., 4./C_AU, 0., 12., 4.)
  • 1: forward ( 1., 90., 14., 4./C_AU, 90., 2., -2.)
  • 2: 90_perp ( 1., 45., 7., 2./C_AU, 45., 7., 1.)
  • 3: 90_par ( 1., 0., 6., 2./C_AU, 0., 6., -2.)
  • 4: integral (4*pi/3., 180., 40., 8*pi/3./C_AU, 180., 40., 0.)

where C_AU is the speed of light in a.u.

Returns:(K_{raman}, x_{1ram}, x_{2ram}, K_{roa}, x_{1roa},
x_{2roa}, x_{3roa})
pyviblib.calc.spectra.stick_raman_roa(scattering, J_all, contr_raman=0, contr_roa=0, contr_weight=True)

Generate Raman/ROA stick spectra.

Parameters:
  • scattering – scattering type (see inv_coeffs_raman_roa())
  • J_all – the five molecular invariants (zero-based ndarray), shape: (NFreq, 5) with NFreq being the number of vibrations returned by RamanROATensors.J_all()
  • contr_raman

    defines how the Raman intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of a2
    • 2 – only the contribution of b2
  • contr_roa

    defines how the ROA intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of aG’
    • 2 – only the contribution of b2G
    • 3 – only the contribution of b2A
  • contr_weight – whether to weight the Raman and ROA intensities. Only relevant if 0 < contr_raman or 0 < contr_roa
Returns:

X, raman_Y, roa_Y, degcirc_Y. Each of these arrays is zero-based and of the length NFreq.

pyviblib.calc.spectra.fit_raman_roa(freqs, J_all, scattering=0, N_G=6, FWHM_is=3.5, FWHM_anis=10.0, FWHM_inst=7.0, startx=None, endx=None, stepx=None, X=None, contr_raman=0, contr_roa=0, contr_weight=True)

Generate Raman/ROA spectra as curves.

Parameters:
  • freqs (one-based ndarray) – wavenumbers of vibrations in ascending order
  • J_all – invariants in approptiate units returned for instance by RamanROATensors.J_all()
  • scattering – scattering type (see inv_coeffs_raman_roa())
  • N_G – number of Gauss functions
  • FWHM_is – full width at half maximum for the isotropic invariants a^2 and b^2
  • FWHM_anis – full width at half maximum for the anisotropic invariants aG', b^2_G, b^2_A
  • FWHM_inst – full width at half maximum for the Gaussian instrument profile
  • startx – start wavenumber of the fitting interval
  • endx – end wavenumbers of the fitting interval if the latter two keyword arguments are None, they are generated automatically
  • stepx – step in cm**(-1), the internal default value is used unless given
  • X – use these x values instead of generating it
  • contr_raman

    defines how the Raman intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of a2
    • 2 – only the contribution of b2
  • contr_roa

    defines how the ROA intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of aG’
    • 2 – only the contribution of b2G
    • 3 – only the contribution of b2A
  • contr_weight – whether to weight the Raman and ROA intensities. Only relevant if 0 < contr_raman or 0 < contr_roa
Returns:

X, raman_Y, roa_Y, degcirc_Y

pyviblib.calc.spectra.make_acp(J, Lx, V)

Generate atomic contribution patterns (ACPs).

Parameters:
  • J – molecular invariant (one-based ndarray), shape: (1 + Natoms, 1 + Natoms), Natoms is the number of atoms
  • Lx – cartesian excursions for the vibration (one-based ndarray), shape: (1 + Natoms, 4)
  • V – correspondent V-tensor (one-based ndarray), shape: (1 + Natoms, 4, 1 + Natoms, 4)
Returns:

the ACPs as a one-based ndarray

pyviblib.calc.spectra.calc_spectra_nmol(mols, representation=0, voa=0, composition=None, startx=None, endx=None, stepx=None, X=None, ab=None, scattering=0, N_G=6, FWHM_is=3.5, FWHM_anis=10.0, FWHM_inst=7.0, fwhm_fit=3.5, groups_to_include=None, contr_raman=0, contr_roa=0, contr_weight=True, distance_terms=0)

Calculate VROA spectra for several molecules.

Parameters:
  • representation – 0 for fitted spectra, 1 for stick spectra
  • voa – 0 for Raman/ROA, 1 for IR/VCD
  • composition – composition of the mixture (zero-based ndarray)
  • startx – start wavenumber of the fitting interval
  • endx – end wavenumbers of the fitting interval if the latter two keyword arguments are None, they are generated automatically
  • stepx – step in cm**(-1), the internal default value is used unless given
  • X – use these x values instead of generating it
  • ab

    if not None, scale the frequencies as follows:

    a = ab[0]
    b = ab[1]
    new_freqs = (a * freqs + b) * freqs
    
  • scattering – scattering type (see inv_coeffs_raman_roa())
  • N_G – number of Gauss functions
  • FWHM_is – full width at half maximum for the isotropic invariants a^2 and b^2
  • FWHM_anis – full width at half maximum for the anisotropic invariants aG', b^2_G, b^2_A
  • FWHM_inst – full width at half maximum for the Gaussian instrument profile
  • fwhm_fit – full width at half maximum of the Lorentz bands
  • groups_to_include – list of indexes of groups whose spectra are to be plotted
  • contr_raman

    defines how the Raman intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of a2
    • 2 – only the contribution of b2
  • contr_roa

    defines how the ROA intensity is to be calculated:

    • 0 – the full intensity
    • 1 – only the contribution of aG’
    • 2 – only the contribution of b2G
    • 3 – only the contribution of b2A
  • contr_weight – whether to weight the Raman and ROA intensities. Only relevant if 0 < contr_raman or 0 < contr_roa
  • distance_terms

    distance terms

    • 0 – usual (local term + non-local term)
    • 1 – local term
    • 2 – non-local term (zero for a2, b2 and aG’)
Returns:

dictionary with the following keys:

  • XY_data
  • mixed_data if composition is given

SourceForge.net Logo