Raman spectroscopy#
The Raman scattering process involves the inelasic scattering of light in the visible range, leaving the molecule in a different vibrational state. The associated scattering amplitude is related to the transition polarizability tensor \(\boldsymbol{\alpha}\), as discussed in more details here. In the case of normal Raman scattering described using the harmonic approximation for molecular vibrations, the transition polarizability associated to a normal mode \(k\) depends on the derivative of the electronic polarizability with respect to the normal coordinate \(Q_k\). As a consequence, only modes along which the polarizability changes will be Raman active.
Spectrum calculation#
Here, we will calculate the Raman spectra of H\(_2\)O, NH\(_3\) and CH\(_4\) and determine how going from O over N to C influences the vibrational frequencies of the Hydrogen stretching modes.
# Import section
import veloxchem as vlx
import py3Dmol as p3d
import sys
from veloxchem.veloxchemlib import bohr_in_angstroms
from matplotlib import pyplot as plt
import numpy as np
basis_set_label = '6-31G'
To compute a meaningful Raman spectrum, the molecular geometries must be optimized at the same level of theory and using the same basis set. Let us assume that we have optimized the geometries at the HF/6-31G theory level and we have these geometries stored in xyz files. Otherwise, we would start by performing the necessary geometry optimizations.
# Prepare molecules and basis sets
# Read optimized geometry from file
h2o_xyz = open("../../data/ir_raman/opt_h2o.xyz", "r").read()
nh3_xyz = open("../../data/ir_raman/opt_nh3.xyz", "r").read()
ch4_xyz = open("../../data/ir_raman/opt_ch4.xyz", "r").read()
h2o = vlx.Molecule.read_xyz_string(h2o_xyz)
h2o_basis = vlx.MolecularBasis.read(h2o, basis_set_label)
nh3 = vlx.Molecule.read_xyz_string(nh3_xyz)
nh3_basis = vlx.MolecularBasis.read(nh3, basis_set_label)
ch4 = vlx.Molecule.read_xyz_string(ch4_xyz)
ch4_basis = vlx.MolecularBasis.read(ch4, basis_set_label)
(a) H2O (b) NH3 (c) CH4
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
jupyter labextension install jupyterlab_3dmol
To compute the Raman spectrum, we need the vibrational frequencies and Raman intensities. The frequencies are eigenvalues of the molecular Hessian, while the Raman intensities are proportional to the gradient of the polarizability. The first step is to determine the Hessian matrix.
# Settings:
scf_settings = {'conv_thresh': 1.0e-5}
method_settings = {}
# Numerical Hessian
hessian_settings = {'do_raman': 'yes', 'print_depolarization_ratio':'yes'}
# We first need to run SCF calculations:
h2o_scf_drv = vlx.ScfRestrictedDriver()
h2o_scf_drv.update_settings(scf_settings, method_settings)
h2o_results = h2o_scf_drv.compute(h2o, h2o_basis)
nh3_scf_drv = vlx.ScfRestrictedDriver()
nh3_scf_drv.update_settings(scf_settings, method_settings)
nh3_results = nh3_scf_drv.compute(nh3, nh3_basis)
ch4_scf_drv = vlx.ScfRestrictedDriver()
ch4_scf_drv.update_settings(scf_settings, method_settings)
ch4_results = ch4_scf_drv.compute(ch4, ch4_basis)
# Create a Hessian driver object, update settings, and compute
h2o_hessian_drv = vlx.ScfHessianDriver()
h2o_hessian_drv.update_settings(method_settings, hessian_settings)
nh3_hessian_drv = vlx.ScfHessianDriver()
nh3_hessian_drv.update_settings(method_settings, hessian_settings)
ch4_hessian_drv = vlx.ScfHessianDriver()
ch4_hessian_drv.update_settings(method_settings, hessian_settings)
# Calculate hessians from scratch
h2o_results = h2o_hessian_drv.compute(h2o, h2o_basis, h2o_scf_drv)
nh3_results = nh3_hessian_drv.compute(nh3, nh3_basis, nh3_scf_drv)
ch4_results = ch4_hessian_drv.compute(ch4, ch4_basis, ch4_scf_drv)
# Or read from checkpoint file:
import h5py
fname = '../../data/ir_raman/raman.h5'
hf = h5py.File(fname, "r")
labels = ['h2o', 'nh3', 'ch4']
i = 0
for driver in [h2o_hessian_drv, nh3_hessian_drv, ch4_hessian_drv]:
label = labels[i]
driver.hessian = np.array(hf.get(label+'_hessian'))
driver.polarizability_gradient = np.array(hf.get(label+'_polgrad'))
driver.dipole_gradient = np.array( hf.get(label+'_dipolegrad'))
i += 1
hf.close()
h2o_hessian_drv.vibrational_analysis(h2o)
Vibrational Analysis
======================
Vibrational Mode 1
----------------------------------------------------
Harmonic frequency: 1736.87 cm**-1
Reduced mass: 1.0917 amu
Force constant: 1.9404 mdyne/A
IR intensity: 123.0358 km/mol
Raman activity: 10.6274 A**4/amu
Parallel Raman: 97.0717 A**4/amu
Perpendicular Raman: 38.4408 A**4/amu
Depolarization ratio: 0.3960
Normal mode:
X Y Z
1 O -0.0554 -0.0410 0.0290
2 H 0.5217 0.0427 -0.4725
3 H 0.3571 0.6079 0.0128
Vibrational Mode 2
----------------------------------------------------
Harmonic frequency: 3988.20 cm**-1
Reduced mass: 1.0372 amu
Force constant: 9.7200 mdyne/A
IR intensity: 2.9544 km/mol
Raman activity: 90.1326 A**4/amu
Parallel Raman: 943.7930 A**4/amu
Perpendicular Raman: 205.5058 A**4/amu
Depolarization ratio: 0.2177
Normal mode:
X Y Z
1 O 0.0327 0.0242 -0.0171
2 H -0.1273 -0.6466 -0.2544
3 H -0.3919 0.2622 0.5260
Vibrational Mode 3
----------------------------------------------------
Harmonic frequency: 4145.10 cm**-1
Reduced mass: 1.0889 amu
Force constant: 11.0231 mdyne/A
IR intensity: 54.2936 km/mol
Raman activity: 40.0871 A**4/amu
Parallel Raman: 292.0906 A**4/amu
Perpendicular Raman: 219.0680 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 O 0.0158 -0.0544 -0.0467
2 H 0.1681 0.6495 0.2172
3 H -0.4197 0.2144 0.5246
nh3_hessian_drv.vibrational_analysis(nh3)
Vibrational Analysis
======================
Vibrational Mode 1
----------------------------------------------------
Harmonic frequency: 597.47 cm**-1
Reduced mass: 1.2003 amu
Force constant: 0.2525 mdyne/A
IR intensity: 644.8631 km/mol
Raman activity: 4.6612 A**4/amu
Parallel Raman: 58.6274 A**4/amu
Perpendicular Raman: 0.8086 A**4/amu
Depolarization ratio: 0.0138
Normal mode:
X Y Z
1 N 0.0794 0.0113 0.0915
2 H -0.4452 -0.0243 -0.3601
3 H -0.3156 0.0200 -0.4779
4 H -0.3420 -0.1531 -0.4336
Vibrational Mode 2
----------------------------------------------------
Harmonic frequency: 1814.46 cm**-1
Reduced mass: 1.0854 amu
Force constant: 2.1054 mdyne/A
IR intensity: 46.9030 km/mol
Raman activity: 6.2677 A**4/amu
Parallel Raman: 45.6689 A**4/amu
Perpendicular Raman: 34.2517 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 N 0.0081 -0.0767 0.0024
2 H 0.1611 0.7186 -0.1654
3 H -0.3086 0.4013 0.3445
4 H 0.0343 -0.0540 -0.2130
Vibrational Mode 3
----------------------------------------------------
Harmonic frequency: 1814.47 cm**-1
Reduced mass: 1.0854 amu
Force constant: 2.1054 mdyne/A
IR intensity: 46.9030 km/mol
Raman activity: 6.2677 A**4/amu
Parallel Raman: 45.6690 A**4/amu
Perpendicular Raman: 34.2517 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 N -0.0579 -0.0045 0.0508
2 H -0.1089 0.2583 -0.1202
3 H 0.3386 -0.4046 -0.0974
4 H 0.5753 0.2093 -0.4884
Vibrational Mode 4
----------------------------------------------------
Harmonic frequency: 3779.47 cm**-1
Reduced mass: 1.0125 amu
Force constant: 8.5217 mdyne/A
IR intensity: 0.5571 km/mol
Raman activity: 109.6648 A**4/amu
Parallel Raman: 1223.8453 A**4/amu
Perpendicular Raman: 174.5125 A**4/amu
Depolarization ratio: 0.1426
Normal mode:
X Y Z
1 N 0.0122 0.0017 0.0141
2 H 0.3677 -0.1623 -0.4143
3 H -0.3412 -0.4046 0.2305
4 H -0.1966 0.5426 -0.0122
Vibrational Mode 5
----------------------------------------------------
Harmonic frequency: 3983.60 cm**-1
Reduced mass: 1.0991 amu
Force constant: 10.2768 mdyne/A
IR intensity: 16.1431 km/mol
Raman activity: 39.5442 A**4/amu
Parallel Raman: 288.1351 A**4/amu
Perpendicular Raman: 216.1014 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 N -0.0518 0.0536 0.0383
2 H 0.4485 -0.2013 -0.5554
3 H 0.0379 0.0729 -0.0208
4 H 0.2330 -0.6159 0.0444
Vibrational Mode 6
----------------------------------------------------
Harmonic frequency: 3983.60 cm**-1
Reduced mass: 1.0991 amu
Force constant: 10.2768 mdyne/A
IR intensity: 16.1430 km/mol
Raman activity: 39.5443 A**4/amu
Parallel Raman: 288.1354 A**4/amu
Perpendicular Raman: 216.1015 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 N 0.0367 0.0639 -0.0398
2 H -0.1983 0.1137 0.2440
3 H -0.4994 -0.5646 0.2943
4 H 0.1873 -0.4371 0.0143
ch4_hessian_drv.vibrational_analysis(ch4)
Vibrational Analysis
======================
Vibrational Mode 1
----------------------------------------------------
Harmonic frequency: 1516.76 cm**-1
Reduced mass: 1.1820 amu
Force constant: 1.6021 mdyne/A
IR intensity: 21.4495 km/mol
Raman activity: 2.9988 A**4/amu
Parallel Raman: 21.8503 A**4/amu
Perpendicular Raman: 16.3878 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 C 0.1208 0.0145 0.0316
2 H -0.3408 0.2769 0.0631
3 H -0.3837 -0.3044 -0.2252
4 H -0.3008 0.0648 -0.3929
5 H -0.4146 -0.2098 0.1782
Vibrational Mode 2
----------------------------------------------------
Harmonic frequency: 1516.77 cm**-1
Reduced mass: 1.1820 amu
Force constant: 1.6021 mdyne/A
IR intensity: 21.4494 km/mol
Raman activity: 2.9987 A**4/amu
Parallel Raman: 21.8500 A**4/amu
Perpendicular Raman: 16.3875 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 C 0.0074 0.1011 -0.0744
2 H 0.0868 -0.3096 0.4673
3 H -0.2187 -0.2678 0.4258
4 H 0.2984 -0.3257 0.0707
5 H -0.2543 -0.3013 -0.0769
Vibrational Mode 3
----------------------------------------------------
Harmonic frequency: 1516.77 cm**-1
Reduced mass: 1.1820 amu
Force constant: 1.6021 mdyne/A
IR intensity: 21.4495 km/mol
Raman activity: 2.9988 A**4/amu
Parallel Raman: 21.8504 A**4/amu
Perpendicular Raman: 16.3878 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 C 0.0340 -0.0734 -0.0963
2 H -0.3867 0.0305 0.2628
3 H 0.1738 -0.0487 0.3378
4 H 0.0638 0.4618 0.2676
5 H -0.2559 0.4307 0.2793
Vibrational Mode 4
----------------------------------------------------
Harmonic frequency: 1708.70 cm**-1
Reduced mass: 1.0080 amu
Force constant: 1.7339 mdyne/A
IR intensity: 0.0000 km/mol
Raman activity: 39.5165 A**4/amu
Parallel Raman: 287.9329 A**4/amu
Perpendicular Raman: 215.9497 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 C 0.0000 0.0000 -0.0000
2 H -0.4146 0.2792 0.0130
3 H -0.3848 -0.2703 -0.1699
4 H 0.3506 -0.1577 0.3197
5 H 0.4487 0.1487 -0.1628
Vibrational Mode 5
----------------------------------------------------
Harmonic frequency: 1708.70 cm**-1
Reduced mass: 1.0080 amu
Force constant: 1.7339 mdyne/A
IR intensity: 0.0000 km/mol
Raman activity: 39.5165 A**4/amu
Parallel Raman: 287.9329 A**4/amu
Perpendicular Raman: 215.9497 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 C -0.0000 -0.0000 -0.0000
2 H 0.1175 0.1950 -0.4452
3 H -0.0637 -0.1956 0.4557
4 H 0.0127 -0.4426 -0.2322
5 H -0.0664 0.4432 0.2217
Vibrational Mode 6
----------------------------------------------------
Harmonic frequency: 3181.65 cm**-1
Reduced mass: 1.0080 amu
Force constant: 6.0118 mdyne/A
IR intensity: 0.0000 km/mol
Raman activity: 145.7008 A**4/amu
Parallel Raman: 1857.8604 A**4/amu
Perpendicular Raman: 0.0000 A**4/amu
Depolarization ratio: 0.0000
Normal mode:
X Y Z
1 C -0.0000 -0.0000 -0.0000
2 H 0.2537 0.3660 0.2273
3 H 0.3128 -0.3724 -0.1161
4 H -0.3562 -0.1710 0.3063
5 H -0.2103 0.1773 -0.4175
Vibrational Mode 7
----------------------------------------------------
Harmonic frequency: 3296.17 cm**-1
Reduced mass: 1.0993 amu
Force constant: 7.0367 mdyne/A
IR intensity: 37.8001 km/mol
Raman activity: 66.7339 A**4/amu
Parallel Raman: 486.2497 A**4/amu
Perpendicular Raman: 364.6872 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 C -0.0340 -0.0610 0.0584
2 H 0.1655 0.2357 0.1716
3 H -0.0733 0.0606 0.0392
4 H 0.5517 0.2531 -0.4667
5 H -0.2385 0.1776 -0.4405
Vibrational Mode 8
----------------------------------------------------
Harmonic frequency: 3296.17 cm**-1
Reduced mass: 1.0993 amu
Force constant: 7.0367 mdyne/A
IR intensity: 37.8000 km/mol
Raman activity: 66.7338 A**4/amu
Parallel Raman: 486.2491 A**4/amu
Perpendicular Raman: 364.6868 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 C -0.0033 -0.0620 -0.0666
2 H 0.3808 0.5344 0.3245
3 H -0.3633 0.4152 0.1170
4 H -0.1206 -0.0737 0.0854
5 H 0.1425 -0.1372 0.2672
Vibrational Mode 9
----------------------------------------------------
Harmonic frequency: 3296.18 cm**-1
Reduced mass: 1.0993 amu
Force constant: 7.0367 mdyne/A
IR intensity: 37.8000 km/mol
Raman activity: 66.7338 A**4/amu
Parallel Raman: 486.2492 A**4/amu
Perpendicular Raman: 364.6869 A**4/amu
Depolarization ratio: 0.7500
Normal mode:
X Y Z
1 C -0.0844 0.0270 -0.0209
2 H 0.1386 0.2390 0.1385
3 H 0.3920 -0.4859 -0.1592
4 H 0.2417 0.1337 -0.2324
5 H 0.2337 -0.2087 0.5026
def add_broadening(list_ex_energy, list_osci_strength, line_profile='Lorentzian', line_param=10, step=10):
x_min = np.amin(list_ex_energy) - 50
x_max = np.amax(list_ex_energy) + 50
x = np.arange(x_min, x_max, step)
y = np.zeros((len(x)))
# go through the frames and calculate the spectrum for each frame
for xp in range(len(x)):
for e, f in zip(list_ex_energy, list_osci_strength):
if line_profile == 'Gaussian':
y[xp] += f * np.exp(-(
(e - x[xp]) / line_param)**2)
elif line_profile == 'Lorentzian':
y[xp] += 0.5 * line_param * f / (np.pi * (
(x[xp] - e)**2 + 0.25 * line_param**2))
return x, y
# plot the Raman spectra
plt.figure(figsize=(6,4))
h2o_x,h2o_y = h2o_hessian_drv.frequencies, h2o_hessian_drv.raman_intensities
nh3_x,nh3_y = nh3_hessian_drv.frequencies, nh3_hessian_drv.raman_intensities
ch4_x,ch4_y = ch4_hessian_drv.frequencies, ch4_hessian_drv.raman_intensities
h2o_xl, h2o_yl = add_broadening(h2o_x, h2o_y, line_profile='Lorentzian', line_param=10, step=10)
nh3_xl, nh3_yl = add_broadening(nh3_x, nh3_y, line_profile='Lorentzian', line_param=10, step=10)
ch4_xl, ch4_yl = add_broadening(ch4_x, ch4_y, line_profile='Lorentzian', line_param=10, step=10)
plt.plot(h2o_xl, h2o_yl, label='H2O')
plt.plot(nh3_xl, nh3_yl, label='NH3')
plt.plot(ch4_xl, ch4_yl, label='CH4')
plt.xlabel('Wavenumber (cm**-1)')
plt.ylabel('Raman activity (A**4/amu)')
plt.title("Calculated Raman sepctra of H2O, NH3, CH4")
plt.legend()
plt.tight_layout(); plt.show()

Let’s compare the H-stretching modes:
# To animate the normal mode we will need both the geometry and the displacements
def get_normal_mode(molecule, normal_mode):
elements = molecule.get_labels()
coords = molecule.get_coordinates() * bohr_in_angstroms() # To transform from au to A
natm = molecule.number_of_atoms()
vib_xyz = "%d\n\n" % natm
nm = normal_mode.reshape(natm, 3)
for i in range(natm):
# add coordinates:
vib_xyz += elements[i] + " %15.7f %15.7f %15.7f " % (coords[i,0], coords[i,1], coords[i,2])
# add displacements:
vib_xyz += "%15.7f %15.7f %15.7f\n" % (nm[i,0], nm[i,1], nm[i,2])
return vib_xyz
h2o_h1 = get_normal_mode(h2o, h2o_hessian_drv.normal_modes[-2])
h2o_h2 = get_normal_mode(h2o, h2o_hessian_drv.normal_modes[-1])
nh3_h1 = get_normal_mode(nh3, nh3_hessian_drv.normal_modes[-3])
nh3_h2 = get_normal_mode(nh3, nh3_hessian_drv.normal_modes[-2])
nh3_h3 = get_normal_mode(nh3, nh3_hessian_drv.normal_modes[-1])
ch4_h1 = get_normal_mode(ch4, ch4_hessian_drv.normal_modes[-4])
ch4_h2 = get_normal_mode(ch4, ch4_hessian_drv.normal_modes[-3])
ch4_h3 = get_normal_mode(ch4, ch4_hessian_drv.normal_modes[-2])
ch4_h4 = get_normal_mode(ch4, ch4_hessian_drv.normal_modes[-1])
print("These are the two H stretching modes of H2O.")
view = p3d.view(viewergrid=(1,2), width=300, height=200, linked=True)
view.addModel(h2o_h1, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,0))
view.addModel(h2o_h2, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,1))
view.setViewStyle({"style": "outline", "width": 0.05})
view.setStyle({"stick":{},"sphere": {"scale":0.25}})
view.animate({'loop': 'backAndForth'})
view.rotate(120, "x")
view.zoomTo()
view.show()
These are the two H stretching modes of H2O.
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
jupyter labextension install jupyterlab_3dmol
print("These are the three H stretching modes of NH3.")
view = p3d.view(viewergrid=(1,3), width=600, height=200, linked=True)
view.addModel(nh3_h1, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,0))
view.addModel(nh3_h2, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,1))
view.addModel(nh3_h3, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,2))
view.setViewStyle({"style": "outline", "width": 0.05})
view.setStyle({"stick":{},"sphere": {"scale":0.25}})
view.animate({'loop': 'backAndForth'})
view.rotate(90, "x")
view.zoomTo()
view.show()
These are the three H stretching modes of NH3.
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
jupyter labextension install jupyterlab_3dmol
print("These are the four H stretching modes of CH4.")
view = p3d.view(viewergrid=(1,4), width=600, height=200, linked=True)
view.addModel(ch4_h1, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,0))
view.addModel(ch4_h2, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,1))
view.addModel(ch4_h3, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,2))
view.addModel(ch4_h4, "xyz", {'vibrate': {'frames':10,'amplitude':0.75}}, viewer=(0,3))
view.setViewStyle({"style": "outline", "width": 0.05})
view.setStyle({"stick":{},"sphere": {"scale":0.25}})
view.animate({'loop': 'backAndForth'})
view.zoomTo()
view.show()
These are the four H stretching modes of CH4.
You appear to be running in JupyterLab (or JavaScript failed to load for some other reason). You need to install the 3dmol extension:
jupyter labextension install jupyterlab_3dmol