Next topic

mmf.plotting

This Page

mmf.physics.seq

ScarfII(*varargin, **kwargs) Scarf II (hyperbolic) potential:
S(*varargin, **kwargs) Scarf II (hyperbolic) potential:
HO_radial(*varargin, **kwargs) Radial Harmonic Oscillator.
R(x, n, a, b) Return the Romanovski polynomial

Inheritance diagram for mmf.physics.seq:

Inheritance diagram of mmf.physics.seq

Numerical Solutions to the Schrodinger equation.

This module contains some exact solutions to the Schrodinger equation for use in testing.

class mmf.physics.seq.ScarfII(*varargin, **kwargs)[source]

Bases: mmf.objects.StateVars

Scarf II (hyperbolic) potential:

ScarfII(A=1,
        B=0,
        a=1,
        m=1,
        hbar=1)

V(x) = \frac{\hbar^2}{2m a^2}\Biggl[A^2 + \left(B^2 - A^2 -
   A\right)
   \sech^{2}\frac{x}{a} + B\left(2A + 1\right)
   \sech\frac{x}{a}\; \tanh \frac{x}{a}\Biggr].

Notes

The depth of the potential is controlled by A while the asymmetry is controlled by B (though the spectrum is unchanged).

from mmf.physics.seq import ScarfII
x = np.linspace(-3,3,100)
for B in [-1,0,1,2]:
   V = ScarfII(A=1, B=B)
   plt.plot(x, V(x), label='B=%i'%(B,))

plt.legend()

(Source code, png, hires.png, pdf)

../../../_images/mmf-physics-seq-1.png

Examples

from mmf.physics.seq import ScarfII
V = ScarfII(A=3, B=2)
x = np.linspace(-10, 10, 200)
v = V(x)
En = V.E()
plt.plot(x, v, 'k')
for n, E in enumerate(En):
   x_ = x[np.where(v <= E)][[0, -1]]
   line = plt.plot(x_, x_*0 + E, lw=2)
   colour = line[0].get_color()
   psi = V.psi(n)(x)
   psi = psi/np.max(abs(psi))
   plt.plot(x, 2*psi, colour + '-')

(Source code, png, hires.png, pdf)

../../../_images/mmf-physics-seq-2.png

Attributes

class mmf.physics.seq.S(*varargin, **kwargs)[source]

Bases: mmf.objects.StateVars

Scarf II (hyperbolic) potential:

S(A=1,
  B=0,
  a=1,
  m=1,
  hbar=1)
ScarfII(A=1,
        B=0,
        a=1,
        m=1,
        hbar=1)

V(x) = \frac{\hbar^2}{2m a^2}\Biggl[A^2 + \left(B^2 - A^2 -
   A\right)
   \sech^{2}\frac{x}{a} + B\left(2A + 1\right)
   \sech\frac{x}{a}\; \tanh \frac{x}{a}\Biggr].

Notes

The depth of the potential is controlled by A while the asymmetry is controlled by B (though the spectrum is unchanged).

from mmf.physics.seq import ScarfII
x = np.linspace(-3,3,100)
for B in [-1,0,1,2]:
   V = ScarfII(A=1, B=B)
   plt.plot(x, V(x), label='B=%i'%(B,))

plt.legend()

(Source code, png, hires.png, pdf)

../../../_images/mmf-physics-seq-3.png

Examples

from mmf.physics.seq import ScarfII
V = ScarfII(A=3, B=2)
x = np.linspace(-10, 10, 200)
v = V(x)
En = V.E()
plt.plot(x, v, 'k')
for n, E in enumerate(En):
   x_ = x[np.where(v <= E)][[0, -1]]
   line = plt.plot(x_, x_*0 + E, lw=2)
   colour = line[0].get_color()
   psi = V.psi(n)(x)
   psi = psi/np.max(abs(psi))
   plt.plot(x, 2*psi, colour + '-')

(Source code, png, hires.png, pdf)

../../../_images/mmf-physics-seq-4.png

Attributes

State Variables:

A

Dimensionless parameter (depth).

B

Dimensionless parameter (controls asymmetry).

a

Length scale

m

<no description>

hbar
class mmf.physics.seq.HO_radial(*varargin, **kwargs)[source]

Bases: mmf.objects.StateVars

Radial Harmonic Oscillator.

HO_radial(d=3,
          hbar=1,
          m=1,
          omega=1)

Notes

In d-dimensions, the Hamiltonian is simply a sum over d independent oscillators, \smash{\op{H}}^{(d)} = \sum_{i=1}^{d}
\op{H}_{i}, where each \op{H}_{i} is an independent oscillator in coordinate x_{i}. The energies are additive and thus:

E_{n_{1},n_{2},\cdots,n_{d}} = \left(\sum_{i} n_{i} +
\frac{d}{2}\right)\hbar\omega.

It is useful to consider this as a central potential problem in the coordinate r. In this case we recover a 1-d problem with an additional set of quantum numbers corresponding to angular momentum l. This follows from the expression of the Laplacian in d-dimensions:

\nabla^{\field{R}^d} = \frac{1}{r^{d-1}}\diff{}{r}r^{d-1}\diff{}{r} +
\frac{1}{r^2}\nabla^{S^{d-1}}

where \nabla^{S^{d-1}} is the Laplacian defined on the unit sphere S^{d-1}. This can be brought into a simple form by introducing the radial momentum operator:

\op{p}_{r} = \frac{-i}{r^{(d-1)/2}}\pdiff{}{r}r^{(d-1)/2}.

In terms of this we have:

\frac{1}{r^{d-1}}\diff{}{r}r^{d-1}\diff{}{r} = -\op{p}_{r}^2 -
\frac{(d-1)(d-3)}{4r^2}.

Finally, introducing the appropriate d-dimensional spherical harmonics \mathcal{Y}_{l}(\vect{r}) we have

\frac{1}{r^2}\nabla^{S^{d-1}}\mathcal{Y}_{l}(\vect{r}) =
-l(l+d-2)\mathcal{Y}_{l}(\vect{r}).

Combining these gives the effective radial Hamiltonianfootnote{One can define raising a lowering operators to analyse this~cite{Arik:2008ix}.}

\begin{aligned}
   \smash{\op{H}}^{(d)} &= \frac{1}{2m}\left(
   \op{p}_{r}^2 + \frac{\nu_{l}^2 - \frac{1}{4}}{\op{r}^2}
   + m^2\omega^2 \op{r}^2 \right), &
   \nu_{l} &= l + \frac{d-2}{2}.
\end{aligned}

The energy levels E_{l,j} and degeneracies g_{l,j} are thus

\begin{aligned}
   E_{l,j} &= \left(l + 2j + \frac{d}{2}\right)\hbar\omega, &
   g_{l,j} &= \binom{d+l+2j-1}{l+2j}
         = \frac{\Gamma(d+l+2j)}{\Gamma(l+2j+1)\Gamma(d)}.
\end{aligned}

with normalized radial eigenfunctions are \Psi^{(l)}_{j}(r) =
R_{j}(r/L) where

R^{l}_{j}(r) = (-1)^{j}\sqrt{\frac{2\Gamma(j+1)}{\Gamma(l+j+d/2)}}
e^{-r^2/2} r^{\nu+1/2}L_{j}^{\nu_l}(r^2)

and L_{j}^{\nu_l} are the generalized Laguerre polynomials.

Attributes

mmf.physics.seq.R(x, n, a, b)[source]

Return the Romanovski polynomial

R_{m}^{(a,b)}(x) &= \frac{1}{w^{(a,b)}(x)}
    \diff[m]{}{x}(1+x^2)^{m}w^{(a,b)}(x)\\
w^{(a,b)}(x) &= (1+x^2)^{-a}e^{b\tan^{-1}(x)}

Notes

Computed using the recurrence

R_{m}^{(a,b)}(x) &= [b + 2(m - a)x]R_{m-1}^{(a,b)}(x) 
   + 2(m - 1)(m - a)(1+x^2)R_{m-2}^{(a-1,b)}(x),\\
R_{0}^{(a,b)} &= 1,\\
R_{1}^{(a,b)} &= 2(1-a)x + b.