Basisbank¶
The braketlab.basisbank
-module provides many standard basis functions used in quantum theory. Here you'll find a brief description of the available functions and examples on how to instantiate them.
The quantum harmonic oscillator¶
The Hamiltonian for the harmonic oscillator is \begin{equation} \hat{H}_0 = \frac{\hat{p}_x^2}{2m_e}+\frac{1}{2}m_e\omega_\text{HO}^2\hat{x}^2, \end{equation} where $\hat{p}_x$ is the $x$ component of the momentum operator, $m_e$ is the electron mass and $\omega_\text{HO}$ is the classical frequency of the HO.
The harmonic oscillator wave functions have the form \begin{equation} \psi_n(x) = N_n H_n(x/x_0) e^{-\frac{1}{2}(x/x_0)^2}, \qquad n=0,1,2,\ldots \end{equation} where $N_n$ is a normalization constant, \begin{equation} N_n = \sqrt{\frac{1}{2^nn!\pi^{1/2}x_0}} \end{equation} and the functions $H_n(x/x_0)$ are Hermite polynomials, the first three of which are given by \begin{align} &H_0(y) = 1 \\ &H_1(y) = 2y \\ &H_2(y) = 4y^2 - 2 \end{align} The Hermite polynomials satisfy the recurrence relation \begin{equation} H_{n+1}(y) = 2yH_n(y) - \frac{\text{d}H_n(y)}{\text{d}y} = 2yH_n(y) - 2nH_{n-1}(y). \end{equation}
These functions can be obtained from basisbank as follows:
import braketlab as bk
import numpy as np
n = 10
p10 = bk.basisbank.get_harmonic_oscillator_function(n)
p10
p10.view(web = True)
The hydrogen atom¶
Disclaimer: there are currently two problems affecting the Hydrogen-orbitals in BraketLab, set to be fixed in a future release. First, they are not properly normalized. Second, they do not visualize on Evince easily, due to how the imaginary part is treated in sympy.
The radial Schrödinger equation for the Hydrogen atom is \begin{equation} \left( -\frac{\hslash^2}{2\mu}\frac{\text{d}^2}{\text{d}r^2} + V_\ell(r) \right) f_{n\ell}(r) = \varepsilon_{n\ell} f_{n\ell}(r), \qquad f_{n\ell}(r) = rR_{n\ell}(r), \end{equation} where \begin{equation} V_\ell(r) = V(r) + \frac{\ell(\ell +1)\hslash^2}{2\mu r^2} \end{equation} This equation can be solved analytically to give the real-valued functions \begin{equation} R_{n\ell}(r) = \sqrt{\left(\frac{2}{na}\right)^3 \frac{(n-\ell-1)!}{2n[(n+\ell)!]^3}} \left( \frac{2}{na} r \right)^\ell {\cal L}^{2\ell+1}_{n+1}(\frac{2}{na}r) e^{-r/na} \end{equation} where \begin{equation} a = \frac{4\pi\epsilon_0\hslash^2}{\mu e^2} \end{equation} has the dimension of length and the associated Laguerre polynomials are defined by \begin{equation} {\cal L}_i^j(x) = \frac{e^xx^{-j}}{i!} \frac{\text{d}^i}{\text{d}x^i} (e^{-x}x^{i+j}) \end{equation} for non-negative integers $i$ and $j$.
The full solution for the hydrogen atom can be cast as a product of the radial function and the spherical harmonics. From basisbank, you obtain these with
n = 1
l = 0
m = 0
s = bk.basisbank.get_hydrogen_function(n,l,m)
s
s.ket_sympy_expression
Solid harmonic Gaussian functions¶
Gaussian-type orbitals are extensively used in quantum chemistry. You may load them from Basisbank with
a = 1.0
l = 4
m = 2
d3 = bk.basisbank.get_gto(a, l, m, position = np.array([0,0,1.0]))
d3.ket_sympy_expression
d3.view()
Slater type orbitals¶
a = 1.0 #exponent
w = 1 #prefactor/weight
n = 1
l = 4
m = 2
d3 = bk.basisbank.get_sto(a,w, n, l, m)
d3.ket_sympy_expression