An Introduction to Rydberg atoms with ARC
In this notebook:
- Configure Notebook
- Rydberg atom energy levels
- Rydberg atom wavefunctions
- Matrix elements for atom interaction with electromagnetic radiation
- Rydberg lifetimes and black-body radiation induced transitions
- Creating highly excited states
- Rydberg atom energy level modifications in DC electric fields (Stark shifts)
Tuning the interactions with external fields (Förster resonance)
- Atomic vapour properties
- Advanced Rydberg examples
- Advanced use of ARC package: interfacing and expansions
Further reading:
For concise general overview of the current research directions in the Rydberg physics, and explanation of elementary underlying ideas, check new eBook below. It was made specifically for final year undergraduates and a new PhD reserchers and uses new electronic format to bring interactive figures for exploring different regimes and building intuition about processes (made with the help of ARC package):
The eBook above is interactive in web version and EPUB version if used with EPUB3 supporting device (e.g. iBooks on Mac/iPad, Gitden Reader in Android, or Document viewer Xreader on Linux). Under Windows we recommend using web version since the best currently available reader in Windows (Calibre) doesn't have full support for EPUB3.
Configure Notebook
This notebook uses ARC package (Alkali Rydberg Calculator) that can be installed from Python Package Index (PyPI) with pip
pip install ARC-Alkali-Rydberg-Calculator
.
For more details see GitHub. This installation can be done from notebook too, by prepending ! to the previous command.
!pip install ARC-Alkali-Rydberg-Calculator
# Configure the matplotlib graphics library and configure it to show
# figures inline in the notebook
%matplotlib inline
import matplotlib.pyplot as plt # Import library for direct plotting functions
import numpy as np # Import Numerical Python
from IPython.core.display import display, HTML #Import HTML for formatting output
# NOTE: Uncomment following lines ONLY if you are not using installation via pip
# import sys, os
# rootDir = '/path/to/arc/directory' # e.g. '/Users/Username/Desktop/ARC-Alkali-Rydberg-Calculator'
# sys.path.insert(0,rootDir)
# import sys, os
# sys.path.insert(0,"..")
import arc
from arc import * #Import ARC (Alkali Rydberg Calculator)
Rydberg Atom Energy Levels
Rydberg states are highly excited states of the outer valence electron where the properties can be scaled in terms of the principal quantum number,
. Originally observed in the spectral lines of hydrogen, the binding energy of the Rydberg series are given by
where
are the quantum numbers, is the Rydberg constant and is the quantum defect. This defect describes the increase in binding energy for an alkali atom with respect to Hydrogen due to penetration and polarisation of the closed inner electron shells, which is most significant for the states with highly elliptical orbits, and can be neglected for states with. The quantum defects are parameterised via
with coefficients
taken from measured energy levels.
The mass-corrected Rydberg constant is given by
where is the electron mass andis the atomic mass of the nucleus, with
To demonstrate the effect of quantum defects, the graph below shows the energy levels of Cesium highlighting the effect of the quantum defect between the different
manifolds.
#Load parameters for Caesium
atom=Caesium()
nmin=6 #Minimum n
nmax=60 #Maximum n
lmin=0 #Minimum l
lmax=3 #Maxmium l
#Plot Energy Levels of Cesium
levels = LevelPlot(atom)
levels.makeLevels(nmin,nmax,lmin,lmax)
levels.drawLevels()
levels.showPlot()
# plot is interactive when called outside the IPython notebook (e.g. from Python program)
n=np.arange(6,100,1)
#Plot Quantum Defects of Cs
fig, axes = plt.subplots(1, 1, figsize=(7,5))
axes.plot(n,atom.getQuantumDefect(n,0,0.5),'.',label="$nS_{1/2}$")
axes.plot(n,atom.getQuantumDefect(n,1,1.5),'.',label="$nP_{3/2}$")
axes.plot(n,atom.getQuantumDefect(n,2,2.5),'.',label="$nD_{5/2}$")
axes.plot(n,atom.getQuantumDefect(n,3,3.5),'.',label="$nF_{7/2}$")
axes.legend(loc=0)
axes.set_xlabel('$n$')
axes.set_ylabel('Quantum Defect $\delta_{n\ell j}$')
axes.set_title('Cs Quantum Defects')
Text(0.5, 1.0, 'Cs Quantum Defects')
Rydberg atom wavefunctions
Rydberg atom wavefunctions can be obtained from Schrödinger's equation,
with reduced mass
anda model potential given by
which at long range gives a
Colomb potential and at short range accounts for the finite size of the core, with radial charge given by
Wavefunctions are calculated by numerical intergration of the model potential using parameters
and taken from Marinescu et.al. Potential includes spin-orbit interaction , whereis fine structure constant.
Rydberg atoms are large, with an average size
, which leads to extremely large dipole matrix elements making them ideal for exploitng strong long range interactions. To illustrate this point, the probabiltiy distribution for atomic wavefunctions are plotted below for increasing which highlights the dramatic scaling with atoms approaching m in size for compared tofor the groundstate!
atom=Caesium()
pqn = [15,30,60] # principal quantum numbers of the states
colors = ["b","r","g"]
l = 0 # S state
j = 0.5 # J = 1/2
plotLegend = []
for i in range(len(pqn)):
n = pqn[i]
step = 0.001
a1,b1 = atom.radialWavefunction(l,0.5,j,\
atom.getEnergy(n, l, j)/27.211,\
atom.alphaC**(1/3.0),\
2.0*n*(n+15.0), step)
legendInfo, = plt.plot(a1,(b1)*(b1),\
"-",lw=2,color = colors[i], label = ("n = %d" % n) )
plotLegend.append(legendInfo)
plt.legend(handles=plotLegend)
plt.xlabel(r"Distance from nucleus $r$ ($a_0$)")
plt.ylabel(r"$\vert rR(r)\vert^2$")
plt.show()
Matrix Elements
Dipole Matrix Elements
Relevant properties of the Rydberg states can be derived through evaluation of dipole matrix elements. Recalling the separability of the atomic wavefunctions into radial and spherical components
in the uncoupled basis, the dipole matrix elements can be expressed as
where
corresponds to andtransitions respectively with the reduced matrix element
where round braces denote Wigner-
symbols and the radial matrix element is evaluated from
For highly excited states the hyperfine-structure splitting becomes negligible, however the fine-structure splitting
means the relevant basis is then fine structure basis () with matrix elements
and the reduced matrix element equal to
where the curly braces denote a Wigner-
symbol.
The cell below demonstrates how to extract the relevent matrix elements for different transitions and highlights the extremely large matrix elements of the Rydberg states for transitions to neighbouring states, which scale
due to the large electron radius shown above. Thus whilst the alkali line transitions have elements of order , the Rydberg states can havemaking them ideally suited for exploiting strong, long range interactions as will be explored below.
Quadrupole Matrix Elements
For accurate interaction potentials at short range, it is necessary to also consider quadrupole matrix elements for the Rydberg atoms. getQuadrupoleMatrixElement
#Cs D2 Transition 6S_{1/2}-->6P_{3/2}
n1=6
l1=0
j1=0.5
mj1=0.5
n2=6
l2=1
j2=1.5
mj2=1.5
q=+1
print("Cs D2 Transition 6S_{1/2}-->6P_{3/2}")
print("====================================")
#Radial Matrix element R_{nlj\rightarrown'l'j'}
print("R_{nlj-n'l'j'} = %.3f ea_0" % atom.getRadialMatrixElement(n1,l1,j1,n2,l2,j2))
#Reduced Matrix Element <l||er||l'>
print("<l||er||l'> = = %.3f ea_0" % atom.getReducedMatrixElementL(n1,l1,j1,n2,l2,j2))
#Reduced Matrix Element <j||er||j'>
print("<j||er||j'> = %.3f ea_0" % atom.getReducedMatrixElementJ(n1,l1,j1,n2,l2,j2))
#Angular Coupling
print("<nljmj|er|n'l'j'mj'> = %.3f ea_0\n" %\
atom.getDipoleMatrixElement(n1,l1,j1,mj1,n2,l2,j2,mj2,q))
#Cs Rydberg Transition 60S_{1/2}-->60P_{3/2}
n1=60
l1=0
j1=0.5
mj1=0.5
n2=60
l2=1
j2=1.5
mj2=1.5
q=+1
print("Cs Rydberg Transition 60S_{1/2}-->60P_{3/2}")
print("===========================================")
#Radial Matrix element R_{nlj\rightarrown'l'j'}
print("R_{nlj-n'l'j'} = %.3f ea_0" % atom.getRadialMatrixElement(n1,l1,j1,n2,l2,j2))
#Reduced Matrix Element <l||er||l'>
print("<l||er||l'> = = %.3f ea_0" % atom.getReducedMatrixElementL(n1,l1,j1,n2,l2,j2))
#Reduced Matrix Element <j||er||j'>
print("<j||er||j'> = %.3f ea_0" % atom.getReducedMatrixElementJ(n1,l1,j1,n2,l2,j2))
#Angular Coupling
print("<nljmj|er|n'l'j'mj'> = %.3f ea_0\n" % \
atom.getDipoleMatrixElement(n1,l1,j1,mj1,n2,l2,j2,mj2,q))
Cs D2 Transition 6S_{1/2}-->6P_{3/2} ==================================== R_{nlj-n'l'j'} = -5.477 ea_0 <l||er||l'> = = 5.477 ea_0 <j||er||j'> = 6.324 ea_0 <nljmj|er|n'l'j'mj'> = 3.162 ea_0 Cs Rydberg Transition 60S_{1/2}-->60P_{3/2} =========================================== R_{nlj-n'l'j'} = 3563.326 ea_0 <l||er||l'> = = -3563.326 ea_0 <j||er||j'> = -4114.575 ea_0 <nljmj|er|n'l'j'mj'> = -2057.287 ea_0
For low-lying states, Numerov integration is usually not accurate enough in estimating dipole matrix element. Package also has a .csv table for each element, that can be manually expanded if needed, with various literature values for dipole matrix elements. These values are either obtained experimentally, or by more advanced theoretical calculations. If value is existing in the literature, package will by default use the literature value with the smallest error estimate.
atom = Caesium()
hasLiteratureValue,dme, info = atom.getLiteratureDME(6,1,1.5,7,0,0.5)
if hasLiteratureValue:
print("%.2f ea_0" % dme)
# additional information about the source
s = "experiment"
if info[0]==1:
s = "theory"
print(" source = %s\n errorEstimate = %.2f\n comment = %s\n reference = %s\n doi = %s" % \
(s,info[1],info[2],info[3],info[4]))
5.61 ea_0 source = theory errorEstimate = 0.01 comment = scaled, table VI reference = Physical Review A, 60 4476 (1999)
https://arc-alkali-rydberg-calculator.readthedocs.io/en/latest/_static/Rydberg_atoms_a_primer.html
No comments:
Post a Comment