pyviblib.calc.spectra
index
/usr/lib/python2.4/site-packages/pyviblib/calc/spectra.py

Module for Raman/ROA as well as IR/VCD spectra generation.
 
The following classes are exported :
    RamanROATensors        -- manipulating the Raman/ROA tensors
    IRVCDTensors           -- manipulating the IR/VCD tensors
    ExpRamanROAData        -- exp. Raman, ROA or Degree of circularity spectra
 
The following functions are exported :
    contract_A_tensor()    -- contract A tensor with unit tensor of Levi-Civita
    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-sections' 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)

 
Modules
       
os
pyviblib

 
Classes
       
pyviblib.util.misc.PropertiesContainer(__builtin__.object)
ExpRamanROAData
IRVCDTensors
RamanROATensors

 
class ExpRamanROAData(pyviblib.util.misc.PropertiesContainer)
    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 resources.STRINGS_EXPSPECTRA_TYPES
    scattering    -- see resources.STRINGS_SCATTERING_TYPES
    wavenumbers   -- wavenumbers
    raman         -- raman
    roa           -- roa
    degcirc       -- degree of circularity
 
 
Method resolution order:
ExpRamanROAData
pyviblib.util.misc.PropertiesContainer
__builtin__.object

Methods defined here:
__init__(self, filename, laser_power, exposure_time, datatype, scattering)
Constructor of the class.
 
Positional arguments :
filename      -- file name
laser_power   -- laser power (mW)
exposure_time -- exposure time (min)
datatype      -- see resources.STRINGS_EXPSPECTRA_TYPES
scattering    -- see resources.STRINGS_SCATTERING_TYPES
__repr__(self)
Can be used to recreate an object with the same value.
__str__(self)
Details.

Data and other attributes inherited from pyviblib.util.misc.PropertiesContainer:
__dict__ = <dictproxy object>
dictionary for instance variables (if defined)
__weakref__ = <attribute '__weakref__' of 'PropertiesContainer' objects>
list of weak references to the object (if defined)

 
class IRVCDTensors(pyviblib.util.misc.PropertiesContainer)
    Class for manipulating the IR/VCD tensors.
 
For details refer to W. Hug and J. Haesler. Int. J. Quant. Chem. 104:695-715,
2005
 
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
 
The following public methods are exported :
    D()                 -- reduced dipole strength in a.u.
    R()                 -- reduced rotational strength in a.u.
    fit_spectra()       -- generate the IR/VCD/g in the curve representation
    integrated_coeffs() -- integrated absorbtion coefficients and g-factor
    norm_sum_vib()      -- normalized sum of VCD over the vibrations
    norm_sum_rot()      -- normalized sum of VCD over the rotations
 
 
Method resolution order:
IRVCDTensors
pyviblib.util.misc.PropertiesContainer
__builtin__.object

Methods defined here:
D(self, p)
Calculate the reduced dipole strength.
 
Positional arguments :
  p -- number of vibration
 
Returned value is expressed in a.u.
R(self, p)
Calculate the reduced rotational strength.
 
Positional arguments :
  p -- number of vibration
 
Returned value is expressed in a.u.
__init__(self, P, M, Lx, freqs)
Constructor of the class.
 
P, Lx and freqs are required to be set.
If M is not set, it will be filled with zeros.    
 
Positional arguments :
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) with Natoms being 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) with Natoms being the number of atoms
Lx    -- cartesian excursions for the vibrations (one-based ndarray)
         shape : (1 + NFreq, 1 + Natoms, 4) with NFreq being the number of
         vibrations
freqs -- wavenumbers of vibrations in ascending order (one-based ndarray)
         shape : (1 + NFreq,)
__repr__(self)
Can be used to recreate an object with the same value.
fit_spectra(self, ngauss=6, fwhm_fit=8.0, fwhm_inst=8.0, startx=None, endx=None)
Generate the IR/VCD spectra as curves.
 
Refer to W. Hug and J. Haesler. Int. J. Quant. Chem. 104:695-715, 2005
 
Keyword arguments :
scattering -- scattering type (see inv_coeffs_raman_roa())
ngauss     -- number of Gauss functions (default 6)
fwhm_fit   -- full width at half maximum for the bandwidthes (default 8.)
fwhm_inst  -- full width at half maximum for the Gaussian instrument profile
              (default 8.0)
startx     -- start wavenumber of the fitting interval (default None)
endx       -- end wavenumbers of the fitting interval (default None)
              if the latter two keyword arguments are None, they are
              generated automatically
 
Return X, eps_Y, deps_Y, g_Y
integrated_coeffs(self)
Get the integrated absorption coeffcients (in SI) and the g-factor.
 
The shape of the bands is considered to be Lorentzian.
 
See Atkins, 6th edition, p. 459.
 
Return a ndarray with a shape of (NFreq, 3).
norm_sum_rot(self, masses)
Calculate the normalized sum of VCD over the rotations.
 
Positional arguments :
masses     -- masses of the atoms in AMU (one-based ndarray)
              shape : (1 + Natoms,)
norm_sum_vib(self, start_p=1, end_p=None)
Calculate the normalized sum of VCD over the vibrations.
 
Keyword arguments :
start_p    -- start vibration (default 1)
end_p      -- end vibration (default None)
              use None to specify the last

Data and other attributes inherited from pyviblib.util.misc.PropertiesContainer:
__dict__ = <dictproxy object>
dictionary for instance variables (if defined)
__weakref__ = <attribute '__weakref__' of 'PropertiesContainer' objects>
list of weak references to the object (if defined)

 
class RamanROATensors(pyviblib.util.misc.PropertiesContainer)
    Class for manipulating the Raman/ROA tensors.
 
For details refer to W. Hug. Chem. Phys., 264(1):53-69, 2001.
 
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)
    lambda_incident -- wavelength of the incident light
 
The following public methods are exported :
    V_a2()          -- V-tensor for the a2 molecular invariant
    V_b2()          -- V-tensor for the b2 molecular invariant
    V_aG()          -- V-tensor for the aG molecular invariant
    V_b2G()         -- V-tensor for the b2G molecular invariant
    V_b2A()         -- V-tensor for the b2A molecular invariant
    a2()            -- matrix of the a2 molecular invariant
    b2()            -- matrix of the b2 molecular invariant
    aG()            -- matrix of the aG molecular invariant
    b2G()           -- matrix of the b2G molecular invariant
    b2A()           -- matrix of the b2A molecular invariant
    J_all()         -- all five molecular invariants for all vibrations
    J_rot()         -- all five reduced molecular invariants for the rotations
    acp()           -- generate atomic contribution patterns (ACPs)
    intensities()   -- calculate the Raman/ROA intensities
    norm_sum_vib()  -- normalized sum of ROA over the vibrations
    norm_sum_rot()  -- normalized sum of ROA over the rotations
 
 
Method resolution order:
RamanROATensors
pyviblib.util.misc.PropertiesContainer
__builtin__.object

Methods defined here:
J_all(self, units='invariant')
Calculate the five molecular invariants for all vibrations.
 
A row is made up of a2, b2, aG, b2G and b2A.
 
Keyword arguments :
units -- units of the invariants (default 'invariant')
         one of ('invariant', 'cross-section', 'A^4/AMU', 'reduced')
 
Return value is a null-based ndarray with the shape of (NFreq, 5) with
NFreq being the number of vibrations.
J_rot(self, Lx_rot)
Calculate the five reduced molecular invariants for the rotations.
 
A row is made up of a2, b2, aG, b2G and b2A.
 
Return value is a null-based ndarray with the shape of (nrot, 5) with
nrot being the number of rotations.
 
Positional arguments :
Lx_rot -- rotations
          shape : (nrot, 1 + Natoms, 4)
 
Lx_rot can be obtained as follows :
 
Lx_rot = mol.Lx_tr_rot[1 + mol.nrot:]
V_a2(self, i=0)
Calculate the V-tensor for the a2 molecular invariant.
 
Keyword arguments :
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.
V_aG(self, i=0)
Calculate the V-tensor for the aG molecular invariant.
 
Keyword arguments :
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.
V_b2(self, i=0)
Calculate the V-tensor for the b2 molecular invariant.
 
Keyword arguments :
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.
V_b2A(self, i=0)
Calculate the V-tensor for the b2A molecular invariant.
 
Keyword arguments :
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.
V_b2G(self, i=0)
Calculate the V-tensor for the b2G molecular invariant.
 
Keyword arguments :
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.
__init__(self, PP, PM, PQ, Lx, freqs, lambda_incident, A=None)
Constructor 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.
 
Positional arguments :
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) with NFreq being 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
 
Keyword arguments :
A               -- gradients of the non-contracted A tensor
                   (one-based ndarray, default None)
                   shape : (1 + Natoms, 4, 4, 4, 4)
__repr__(self)
Can be used to recreate an object with the same value.
a2(self, p, i=0, units='invariant')
Calculate the a2 molecular invariant.
 
Positional arguments :
p     -- number of vibration
 
Keyword arguments :
i     -- part of the decomposed invariant (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 (default 'invariant')
         one of ('invariant', 'cross-section', 'A^4/AMU')
aG(self, p, i=0, units='invariant')
Calculate the aG molecular invariant.
 
Positional arguments :
p     -- number of vibration
 
Keyword arguments :
i     -- part of the decomposed invariant (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 (default 'invariant')
         one of ('invariant', 'cross-section', 'A^4/AMU')
acp(self, item, p)
Generate atomic contribution patterns (ACPs).
 
See also make_acp().
 
Positional arguments :
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 units of the cross-sections.
b2(self, p, i=0, units='invariant')
Calculate the b2 molecular invariant.
 
Positional arguments :
p     -- number of vibration
 
Keyword arguments :
i     -- part of the decomposed invariant (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 (default 'invariant')
         one of ('invariant', 'cross-section', 'A^4/AMU')
b2A(self, p, i=0, units='invariant')
Calculate the b2A molecular invariant.
 
Positional arguments :
p     -- number of vibration
 
Keyword arguments :
i     -- part of the decomposed invariant (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 (default 'invariant')
         one of ('invariant', 'cross-section', 'A^4/AMU')
b2G(self, p, i=0, units='invariant')
Calculate the b2G molecular invariant.
 
Positional arguments :
p     -- number of vibration
 
Keyword arguments :
i     -- part of the decomposed invariant (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 (default 'invariant')
         one of ('invariant', 'cross-section', 'A^4/AMU')
norm_sum_rot(self, masses, scattering=0)
Calculate the normalized sum of ROA over the rotations.
 
See eq 33 in W. Hug and J. Haesler. Int. J. Quant. Chem. 104:695-715, 2005.
 
Positional arguments :
masses     -- masses of the atoms in AMU (one-based ndarray)
              shape : (1 + Natoms,)
 
Keyword arguments :
scattering -- scattering type (see inv_coeffs_raman_roa()) (default 0)
norm_sum_vib(self, scattering=0, start_p=1, end_p=None)
Calculate the normalized sum of ROA over the 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 W. Hug and J. Haesler. Int. J. Quant. Chem. 104:695-715, 2005.
 
Keyword arguments :
scattering -- scattering type (see inv_coeffs_raman_roa()) (default 0)
start_p    -- start vibration (default 1)
end_p      -- end vibration (default None)
              use None to specify the last

Static methods defined here:
intensities(invars, scattering)
Calculate the Raman/ROA intensities.
 
Positional arguments :
invars     -- array with the five invariants (ndarray)
              shape : (5,)
scattering -- scattering type (see inv_coeffs_raman_roa())

Data and other attributes inherited from pyviblib.util.misc.PropertiesContainer:
__dict__ = <dictproxy object>
dictionary for instance variables (if defined)
__weakref__ = <attribute '__weakref__' of 'PropertiesContainer' objects>
list of weak references to the object (if defined)

 
Functions
       
Kp(nu_0, nu)
Factor appearing in the calculating of the Raman/ROA cross-sections.
 
Kp = 10^7 / 9 * pi^2 * MU_VACUUM^2 * C^4 * (nu_0 - nu)^3 * nu_0
 
For details refer to W. Hug. Chem. Phys., 264(1):53-69, 2001.
 
Positional arguments :
nu_0  -- wavenumber of the incident light in cm**(-1)
nu    -- shift in cm**(-1)
 
Return value is expressed in SI units.
contract_A_tensor(A)
Contract A tensor with antisymmetric unit tensor of Levi-Civita.
 
PQ[m, n] = sum_{l, k} eps[m, l, k] * A[l, k, n]
 
Positional arguments :
A -- A tensor (one-based ndarray)
     shape : (1 + Natoms, 4, 4, 4, 4) with Natoms being the number of atoms
 
Return the contracted tensor PQ with a shape of (1 + Natoms, 4, 4, 4).
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.
.
  Positional arguments :
  J -- molecular invariant
       one of ('a2', 'b2', 'aG', 'b2G', 'b2A')
conv_cross_sections(lambda_incident, nu)
Conversion factor from J_ab in a.u. to the cross-sections in A^2 / sr.
 
For details refer to W. Hug. Chem. Phys., 264(1):53-69, 2001.
 
Positional arguments :
lambda_incident -- wavelength of the incident light in nm
nu              -- wavenumber in cm**(-1)
conv_invariant(nu)
Conversion factor from J_ab in a.u. to molecular invariants in SI units.
 
For details refer to W. Hug. Chem. Phys., 264(1):53-69, 2001.
 
Positional arguments :
nu  -- wavenumber in cm**(-1)
conv_units_raman_roa(units, lambda_incident, nu, J)
Conversion factor for a Raman/ROA quantity.
 
Positional arguments :
units           -- units of the molecular invariant
                   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')
f_Qp_i(nu)
Transition integral between two molecular vibrational states |i> and <f|.
 
<f|Qp|i> = sqrt(HBAR / 400 * pi * C * nu)
 
It is assumed that the f <- i is a fundamental transition and that the
vibration is harmonic, i.e. described by the normal mode p. For details
refer to W. Hug. Chem. Phys., 264(1):53-69, 2001.
 
Positional arguments :
nu -- wavenumber in cm**(-1)
 
Return value is expressed in SI units.
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)
Generate Raman/ROA spectra as curves.
 
For details refer to W. Hug and J. Haesler. Int. J. Quant. Chem. 104:695-715,
2005
 
Positional arguments :
freqs      -- wavenumbers of vibrations in ascending order (one-based array)
J_all      -- invariants in approptiate units
              returned by RamanROATensors.J_all()
 
Keyword arguments :
scattering -- scattering type (see inv_coeffs_raman_roa())
N_G        -- number of Gauss functions (default 6)
FWHM_is    -- full width at half maximum for the isotropic invariants
              a2 and b2 (default 3.5)
FWHM_anis  -- full width at half maximum for the anisotropic invariants
              b2, b2G, b2A (default 10.)
FWHM_instr -- full width at half maximum for the Gaussian instrument profile
              (default 7.0)
startx     -- start wavenumber of the fitting interval (default None)
endx       -- end wavenumbers of the fitting interval (default None)
              if the latter two keyword arguments are None, they are generated
              automatically
 
Return X, raman_Y, roa_Y, degcirc_Y.
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.
 
RAMAN = K_raman * (x1_ram * a2 + x2_ram * b2)
ROA   = K_roa   * (x1_roa * aG + x2_roa * b2G + x3_roa * b2A)
 
For details refer to W. Hug. Chem. Phys., 264(1):53-69, 2001.
 
Keyword arguments :
scattering -- scattering type (default 0)
 
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 : ICP_perp  (     1.,  45.,  7.,      2./C_AU,  45.,  7.,  1.)
    3 : ICP_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.
 
Return (K_raman, x1_ram, x2_ram, K_roa, x1_roa, x2_roa, x3_roa)
make_acp(J, Lx, V)
Generate atomic contribution patterns (ACPs).
 
For details refer to W. Hug. Chem. Phys., 264(1):53-69, 2001.
 
Positional arguments :
J  -- molecular invariant (one-based ndarray)
      shape : (1 + Natoms, 1 + Natoms) with Natoms being 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)
 
Return the ACPs as a one-based ndarray.
stick_raman_roa(scattering, J_all)
Generate Raman/ROA stick spectra.
 
Positional arguments :
scattering -- scattering type (see inv_coeffs_raman_roa())
J_all      -- the five molecular invariants (null-based ndarray)
              shape : (NFreq, 5) with NFreq being the number of vibrations
              returned by RamanROATensors.J_all()
 
Return X, raman_Y, roa_Y, degcirc_Y. Each of these arrays is null-based and
of the length NFreq.

 
Data
        __all__ = ['RamanROATensors', 'IRVCDTensors', 'ExpRamanROAData', 'contract_A_tensor', 'f_Qp_i', 'Kp', 'conv_units_raman_roa', 'conv_invariant', 'conv_cross_sections', 'conv_A4_AMU', 'inv_coeffs_raman_roa', 'stick_raman_roa', 'fit_raman_roa', 'make_acp']
__author__ = 'Maxim Fedorovsky'

 
Author
        Maxim Fedorovsky