Table Of Contents

This Page

mmf.physics.potentials.two_body

DoubleExponential([v0, r0, mu0, f_v, f_r, ...]) Double Exponential potential:
SoonYongPotential2([v0, A0, V1, A1, hbar, ...]) Extended Potential used by Soon Yong.
Potential([abs_tol, rel_tol, Rmin]) Instances are callable and return compact potentials V(r).
PoschlTeller([v0, r0, hbar, M, abs_tol, ...]) Extended Poschl-Teller (Rosen-Morse) potential used by Soon
BargmannPotential([a, re, hbar, M]) Bargmann Potential where the effective range expansion is exact
PotentialV0R0([abs_tol, rel_tol, Rmin]) Mixing for two-parameter potentials with a parameter v_0 governing the overall strength and a parameter r_0 governing the effective range.
DoubleGaussian([v0, r0, mu0, f_v, f_r, ...]) Double Gaussian potential:
SoonYongPotential Extended Poschl-Teller (Rosen-Morse) potential used by Soon

Inheritance diagram for mmf.physics.potentials.two_body:

Inheritance diagram of mmf.physics.potentials.two_body

Code for extracting properties of two-body systems.

Phase Shifts

We consider here short-range potentials such that for V(r>R) = 0. In this case, the potential may be described in terms of the phase shift of the scattering states.

class mmf.physics.potentials.two_body.DoubleExponential(v0=1.1909151179666313, r0=0.6666666666666666, mu0=None, f_v=-2.0, f_r=2.0, hbar=1.0, M=0.5, abs_tol=1e-12, rel_tol=1e-12, Rmin=0.1)[source]

Bases: mmf.physics.potentials.two_body.PotentialV0R0

Double Exponential potential:

V(r) = -v_0 \frac{\hbar^2}{M} \mu_0^2 \left(
  \exp[-\mu_0 r] + f_v \exp[-\mu_0 r f_r]
\right)

Note that \mu_0 = 4/r_0 here.

>>> v = DoubleExponential(mu0=14.0)
>>> v.get_ar()              # One of Stefano's values to check.
(5852295157..., 0.3099...)
>>> v.set_ainv_re(ainv=0.0, re=1.0)
>>> v.v0, v.r0
(1.190915117966..., 0.9218838916...)

Methods

__call__(r)
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V’(r)/V(r). Should be overloaded for speed and/or
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
get_mu0()
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.
set_ainv_re(ainv, re) Uses a non-linear root-finder to set the parameters v0 and r0 so that the scattering length is 1.0/ainv and the effective range is re.
set_mu0(mu0)
__init__(v0=1.1909151179666313, r0=0.6666666666666666, mu0=None, f_v=-2.0, f_r=2.0, hbar=1.0, M=0.5, abs_tol=1e-12, rel_tol=1e-12, Rmin=0.1)[source]
get_mu0()[source]
mu0
set_mu0(mu0)[source]
class mmf.physics.potentials.two_body.SoonYongPotential2(v0=-1.0, A0=6.2527, V1=0.0, A1=3.0, hbar=0.15915494309189535, M=0.5, kF=1.0, abs_tol=1e-12, rel_tol=1e-12, Rmin=0.1)[source]

Bases: mmf.physics.potentials.two_body.Potential

Extended Potential used by Soon Yong.

V(r) = hbar**2/M*(v0*mu0**2/(cosh(mu0*r))**2
  • V1*mu1**2/(cosh(mu1*r))**2)

His units are set by kF = (3*pi*pi*N)**(1./3.) / L for a box of length L mu = A*kF

As a test, at unitarity, v0 = 1.0 and re = 2/mu for two-body scattering where m here is the individual mass of the particles.

>>> V = SoonYongPotential2(v0=-1.0,V1=0.0)
>>> (a,re) = V.get_ar()
>>> mu = V.A0*V.kF
>>> abs(a) > 6e8
True
>>> abs(re - 2.0/mu) < 1e-9
True

# kF_SY = (3.*pi)**(2./3.) # EFG_SY = kF_SY**2/2.0*3.0/5.0 # EFG = EFG_SY/4.0/pi**2 # re = 2.0*L/(6.2527*(3.0*pi)**(2./3.)) # # V = -2*hbar^2*mu^2*v0/m/cosh(mu*r)**2 # v0 = 1.0 => ainv = 0, re = 2.0/mu # v0 = 1.30119 => kfa_SY = 1.0

Methods

__call__(r)
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V’(r)/V(r).
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.
set_ainv([ainv]) Set the coefficient v0 so that the scattering length is a.
set_ainv_re([ainv, re]) Set the coefficients v0 and A0 so that the scattering length is
set_ainv_re_VV([ainv, re]) Set the coefficients v0 and V1 so that the scattering length is
__init__(v0=-1.0, A0=6.2527, V1=0.0, A1=3.0, hbar=0.15915494309189535, M=0.5, kF=1.0, abs_tol=1e-12, rel_tol=1e-12, Rmin=0.1)[source]
dVdV(r)[source]

Return V’(r)/V(r).

V:=v0/cosh(mu0*r)^2+V1/cosh(mu1*r)^2; with(codegen); C(simplify(diff(V,r)/V),optimized);

set_ainv(ainv=0.0)[source]

Set the coefficient v0 so that the scattering length is a.

>>> V = SoonYongPotential2()
>>> V.v0 = -0.5
>>> V.set_ainv(0.0)
>>> np.allclose(V.v0, -1.0)
True
set_ainv_re(ainv=0.0, re=0.07)[source]

Set the coefficients v0 and A0 so that the scattering length is 1.0/ainv and the effective range is re.

>>> V = SoonYongPotential2(kF=1.0)
>>> V.V1 = 0.0
>>> V.set_ainv_re(ainv=0.0,re=1.0)
>>> np.allclose(V.A0, 2./V.kF)
True
>>> np.allclose(V.v0, -1.0)
True
set_ainv_re_VV(ainv=0.0, re=0.07)[source]

Set the coefficients v0 and V1 so that the scattering length is 1.0/ainv and the effective range is re.

>>> V = SoonYongPotential2()
>>> V.A0 = 6.2527
>>> V.A1 = 3.0
>>> V.v0 = -1.5
>>> V.V1 = 0.6
>>> V.set_ainv_re_VV(ainv=0.0,re=0.0)
>>> np.allclose(V.v0, -1.4789900537007772)
True
>>> np.allclose(V.V1, 0.61932931745837538)
True
class mmf.physics.potentials.two_body.Potential(abs_tol=1e-12, rel_tol=1e-12, Rmin=None)[source]

Bases: object

Instances are callable and return compact potentials V(r).

The class allows one to specify parameters. The important property is that when parameters are modified, the class ensures that the radius returned by R = get_R() is such that V(R) < abs_tol.

The parameter M should be the reduced mass of the system.

Methods

__call__(r) Return V(r).
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V’(r)/V(r). Should be overloaded for speed and/or
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.

Construct the potential object.

Methods

__call__(r) Return V(r).
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V’(r)/V(r). Should be overloaded for speed and/or
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.
__init__(abs_tol=1e-12, rel_tol=1e-12, Rmin=None)[source]

Construct the potential object.

R

Return a radius for which self(R) <= abs_tol.

compute_R()[source]

Return R >= self.Rmin such that abs(self(R)) < self.abs_tol.

dVdV(r)[source]

Return V’(r)/V(r). Should be overloaded for speed and/or accuracy.

gE(E, R=None, method='odeint')[source]

Return k*cot(delta) as a function of energy E for the specified scattering problem by matching at a specified radius R.

get_ar(R=None, method='odeint')[source]

Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range. Integrates the wavefunction to radius R. Answer is independent of R if V(r>R) = 0.

get_deq(f, jac, U0, method='odeint')[source]

Return a “deq” object capable of integrating the system defined by U’ = f(r,U) and diff(U’,U) = jac(r,U) with initial condition U(0.0) = U0.

get_deq_E(E, method='odeint')[source]

Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.

Uses the following system: du’(r) = 2*m/hbar**2*(V(r)-E)*u(r) u’(r) = du(r)

with initial conditions

u(0) = 0.0 du(0) = 1.0

Call (u,du) = deq.integrate(R) to integrate to R. Note that deq maintains state, so subsequent calls must have increasing R.

This is to be used by calling the method integrate(R).

get_deq_ar(method='odeint')[source]

Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.

Uses the following system: u’(r) = du(r) du’(r) = 2*m/hbar**2*V(r)*u(r) uu’(r) = duu(r) duu’(r) = 2*m/hbar**2*(V(r)*uu(r)-u(r))

with initial conditions

u(0) = 0.0 du(0) = 1.0 uu(0) = 0.0 duu(0) = 0.0

Call (u,du,uu,duu) = deq.integrate(R) to integrate to R. Note that deq maintains state, so subsequent calls must have increasing R.

get_deq_gE(E, method='odeint')[source]

Return solution to energy dependent differential equation for (Ap, Am) for finding the s-wave phase shift which is

g(E) = \kappa \frac{A_+ - A_-}{A_+ + A_-}

where E = -\kappa^2 and the s-wave radial wavefunction has the asymptotic form

u(r>R) = A_+ e^{\kappa r} + A_-e^{-\kappa r}.

We solve the radial Schrödinger equation

u''(r) = [\tilde{V}(r) - \tilde{E}]u(r)

where \tilde{V} = 2mV(r)/\hbar^2 and \tilde{E} = 2mE(r)/\hbar^2 with initial conditions u(0)=0, and u'(0)=1 using the following transformation

A_{\pm}(r) = \frac{1}{2}\left(u(r) \pm \frac{u'(r)}{\kappa}\right)
             e^{\mp\kappa r}

These should be constant at the asymptotic values once r>R is larger than the support of the potential. These functions satisfy the following set of coupled first-order equations:

\diff{}{r}
\begin{pmatrix}
  A_+(r)\\
  A_-(r)
\end{pmatrix} = \frac{\tilde{V}(r)}{2\kappa}\begin{pmatrix}
  1 & e^{-2\kappa r}\\
  -e^{2\kappa r} & -1
\end{pmatrix}
\begin{pmatrix}
  A_+(r)\\
  A_-(r)
\end{pmatrix}

with initial conditions A_{\pm}(0) = \pm\frac{1}{2\kappa}.

get_deq_gE_old(E, method='odeint', version=1)[source]

Old version in terms of fp and fm.

version=0

fp' = dfp'
dfp' = (kappa**2 + VV)*fp - V'/V*(kappa*fp - dfp)
fm' = dfm
dfm' = (kappa**2 + VV)*fm + V'/V*(kappa*fm + dfm)

fp(0) = 1/2/kappa
dfp(0) = 1/2
fm(0) = -1/2/kappa
dfm(0) = 1/2

version=1

fp' =  ((2*kappa**2 + VV)*fp + VV*fm)/2/kappa
fm' = -((2*kappa**2 + VV)*fm + VV*fp)/2/kappa

fp(0) = 1/2/kappa
fm(0) = -1/2/kappa

Call (fp,dfp,fm,dfm) = deq.integrate(R) to integrate to R. Note that deq maintains state, so subsequent calls must have increasing R.

This is to be used by calling the method integrate(R).

Examples

>>> V = DoubleGaussian(f_v=0.0)
>>> V._gE_neg(E=-1.0, version=0)
-1.874095805...
>>> V._gE_neg(E=-1.0, version=1)
-1.874095805...
>>> V._gE_neg(E=-1.0, version=2)
-1.874095805...
plot(Rmax=None)[source]

Plot potential.

plot_gE(Erange, Rmax=None)[source]

Plot potential.

plot_u(E, Rmax=None)[source]

Plot radial wavefunction.

class mmf.physics.potentials.two_body.PoschlTeller(v0=-1.0, r0=0.31986181969389227, hbar=1.0, M=0.5, abs_tol=1e-12, rel_tol=1e-12, Rmin=0.1)[source]

Bases: mmf.physics.potentials.two_body.PotentialV0R0

Extended Poschl-Teller (Rosen-Morse) potential used by Soon Yong, Alex, Joe, and Stefano:

V(r) &= v_0 \frac{\hbar^2}{M}\frac{\mu_0^2}{\cosh^2(\mu_0 r)}\\
V(r) &= v_0 \frac{\hbar^2}{M}\frac{4}{r_0^2}\sech^2\frac{2r}{r_0}\\

Units are set by k_F = (3\pi^2N_+)^{1/3}/L for a box of length L and \mu = A k_F. As a test, at unitarity, v_0 = 1 and r_e = 2/\mu_0 for two-body scattering where m here is the individual mass of the particles (M=0.5 by default).

>>> V = PoschlTeller(v0=-1.0)
>>> (a,re) = V.get_ar()
>>> mu = 2.0/V.r0
>>> abs(a) > 6e8
True
>>> abs(re - 2.0/mu) < 1e-9
True

It is resonant (a=\infty) for v_0 = -n(2n - 1) for n \in [1,2,3,...].

>>> def infty(n):
...    return abs(PoschlTeller.a_(-n*(2.0*n-1.0)))
>>> infty(1) > 1e8
True
>>> infty(2) > 1e8
True

At these resonances one has r_e/r_0 = H(2n-1) where H(n) is the Harmonic Number H(n) = \sum_{j=1}^{n}j^{-1}.

>>> def H(n):
...     return (1./numpy.arange(1,n+1)).sum()
>>> def zero(n):
...    return abs(PoschlTeller.re_(-n*(2.0*n-1.0))-H(2*n-1))
>>> zero(1) < 1e-10
True
>>> zero(2) < 1e-10
True

The scattering length is zero at the following points where the effective range diverges. I am not sure yet what these correspond to:

v0 = [0.0,-4.4438633257625231,-12.913326147200248,
-25.421577034814653,-41.955714059176017,
-62.508749555456298,...]

# kF_SY = (3.*pi)**(2./3.)
# EFG_SY = kF_SY**2/2.0*3.0/5.0
# EFG = EFG_SY/4.0/pi**2
# re = 2.0*L/(6.2527*(3.0*pi)**(2./3.))
#
# V = -2*hbar^2*mu^2*v0/m/cosh(mu*r)**2
# v0 = 1.0 => ainv = 0, re = 2.0/mu
# v0 = 1.30119 => kfa_SY = 1.0

Methods

__call__(r)
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V'(r)/V(r)
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
k_cot_delta_analytic(E[, l]) This potential has some analytic properties (see Alex’s thesis and
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.
set_ainv([ainv]) Set the coefficient v0 so that the inverse scattering length is ainv.
set_ainv_re([ainv, re]) Set the coefficients v0 and r0 so that the inverse scattering length is ainv and the effective range is re.
__init__(v0=-1.0, r0=0.31986181969389227, hbar=1.0, M=0.5, abs_tol=1e-12, rel_tol=1e-12, Rmin=0.1)[source]
a_ = <numpy.lib.function_base.vectorize object at 0x9975910>[source]
are_ = <numpy.lib.function_base.vectorize object at 0x99758b0>[source]
dVdV(r)[source]

Return V'(r)/V(r)

V:=v0/cosh(mu0*r)^2+V1/cosh(mu1*r)^2;
with(codegen);
C(simplify(diff(V,r)/V),optimized);
k_cot_delta_analytic(E, l=None)[source]

This potential has some analytic properties (see Alex’s thesis and references therein):

System Message: WARNING/2 (k\cot\delta_k &= -k\frac{\Im{Z}}{\Re{Z}},\\ Z &= \frac{2^{-i\tilde{k}}\Gamma(i\tilde{k})} {\Gamma\left(\frac{1+\lambda + i\tilde{k}}{2}\right) \Gamma\left(\frac{2-\lambda + i\tilde{k}}{2}\right)},\\ \tilde{k} &= k/\mu, \qquad l = \frac{1 + \sqrt{1 - 8v_0}}{2}} )

latex exited with error: [stderr] [stdout] This is pdfTeX, Version 3.1415926-2.3-1.40.12 (TeX Live 2011) restricted \write18 enabled. entering extended mode (./math.tex LaTeX2e <2011/06/27> Babel <v3.8m> and hyphenation patterns for english, dumylang, nohyphenation, ge rman-x-2011-07-01, ngerman-x-2011-07-01, afrikaans, ancientgreek, ibycus, arabi c, armenian, basque, bulgarian, catalan, pinyin, coptic, croatian, czech, danis h, dutch, ukenglish, usenglishmax, esperanto, estonian, ethiopic, farsi, finnis h, french, friulan, galician, german, ngerman, swissgerman, monogreek, greek, h ungarian, icelandic, assamese, bengali, gujarati, hindi, kannada, malayalam, ma rathi, oriya, panjabi, tamil, telugu, indonesian, interlingua, irish, italian, kurmanji, lao, latin, latvian, lithuanian, mongolian, mongolianlmc, bokmal, nyn orsk, polish, portuguese, romanian, romansh, russian, sanskrit, serbian, serbia nc, slovak, slovenian, spanish, swedish, turkish, turkmen, ukrainian, uppersorb ian, welsh, loaded. (/usr/local/texlive/2011/texmf-dist/tex/latex/base/article.cls Document Class: article 2007/10/19 v1.4h Standard LaTeX document class (/usr/local/texlive/2011/texmf-dist/tex/latex/base/size12.clo)) (/usr/local/texlive/2011/texmf-dist/tex/latex/base/inputenc.sty (/usr/local/texlive/2011/texmf-dist/tex/latex/ucs/utf8x.def)) (/usr/local/texlive/2011/texmf-dist/tex/latex/ucs/ucs.sty (/usr/local/texlive/2011/texmf-dist/tex/latex/ucs/data/uni-global.def)) (/usr/local/texlive/2011/texmf-dist/tex/latex/amsmath/amsmath.sty For additional information on amsmath, use the `?’ option. (/usr/local/texlive/2011/texmf-dist/tex/latex/amsmath/amstext.sty (/usr/local/texlive/2011/texmf-dist/tex/latex/amsmath/amsgen.sty)) (/usr/local/texlive/2011/texmf-dist/tex/latex/amsmath/amsbsy.sty) (/usr/local/texlive/2011/texmf-dist/tex/latex/amsmath/amsopn.sty)) (/usr/local/texlive/2011/texmf-dist/tex/latex/amscls/amsthm.sty) (/usr/local/texlive/2011/texmf-dist/tex/latex/amsfonts/amssymb.sty (/usr/local/texlive/2011/texmf-dist/tex/latex/amsfonts/amsfonts.sty)) (/usr/local/texlive/2011/texmf-dist/tex/latex/tools/bm.sty) (/Users/mforbes/Library/texmf/tex/latex/mmf/mmfmath.sty (/usr/local/texlive/2011/texmf-dist/tex/latex/mh/mathtools.sty (/usr/local/texlive/2011/texmf-dist/tex/latex/graphics/keyval.sty) (/usr/local/texlive/2011/texmf-dist/tex/latex/tools/calc.sty) (/usr/local/texlive/2011/texmf-dist/tex/latex/mh/mhsetup.sty)) (/usr/local/texlive/2011/texmf-dist/tex/latex/base/ifthen.sty) (/usr/local/texlive/2011/texmf-dist/tex/latex/doublestroke/dsfont.sty) (/usr/local/texlive/2011/texmf-dist/tex/generic/mathdots/mathdots.sty)) (/usr/local/texlive/2011/texmf-dist/tex/latex/preview/preview.sty) (./math.aux) (/usr/local/texlive/2011/texmf-dist/tex/latex/ucs/ucsencs.def) (/usr/local/texlive/2011/texmf-dist/tex/latex/graphics/graphicx.sty (/usr/local/texlive/2011/texmf-dist/tex/latex/graphics/graphics.sty (/usr/local/texlive/2011/texmf-dist/tex/latex/graphics/trig.sty) (/usr/local/texlive/2011/texmf-dist/tex/latex/latexconfig/graphics.cfg) (/usr/local/texlive/2011/texmf-dist/tex/latex/graphics/dvips.def))) Preview: Fontsize 12pt (/usr/local/texlive/2011/texmf-dist/tex/latex/amsfonts/umsa.fd) (/usr/local/texlive/2011/texmf-dist/tex/latex/amsfonts/umsb.fd) ! Argument of \gather has an extra }. <inserted text> \par l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}} \end{split}\notag\\\begin{split}\end{... Runaway argument? \begin {split}k\cot \delta _k &= -k\frac {\Im {Z}}{\Re {Z}},\\ Z &= \ETC. ! Paragraph ended before \gather was complete. <to be read again> \par l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}} \end{split}\notag\\\begin{split}\end{... ! Missing $ inserted. <inserted text> $ l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}} \end{split}\notag\\\begin{split}\end{... ! Missing \endgroup inserted. <inserted text> \endgroup l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}} \end{split}\notag\\\begin{split}\end{... ! Display math should end with $$. <to be read again> \par l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}} \end{split}\notag\\\begin{split}\end{... ! Extra }, or forgotten \endgroup. <recently read> } l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}} \end{split}\notag\\\begin{split}\end{... ! Misplaced \crcr. \endsplit ->\crcr \egroup \egroup \iftagsleft@ \@xp \lendsplit@ \else \@xp \... l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}}\end{split} \notag\\\begin{split}\end{... ! Extra }, or forgotten \endgroup. \endsplit ->\crcr \egroup \egroup \iftagsleft@ \@xp \lendsplit@ \else \@xp \... l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}}\end{split} \notag\\\begin{split}\end{... ! Extra }, or forgotten \endgroup. \endsplit ->\crcr \egroup \egroup \iftagsleft@ \@xp \lendsplit@ \else \@xp \... l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}}\end{split} \notag\\\begin{split}\end{... ! LaTeX Error: \begin{gather} on input line 16 ended by \end{split}. See the LaTeX manual or LaTeX Companion for explanation. Type H <return> for immediate help. ... l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}}\end{split} \notag\\\begin{split}\end{... ! Missing $ inserted. <inserted text> $ l.22 l = \frac{1 + \sqrt{1 - 8v_0}}{2}}\end{split} \notag\\\begin{split}\end{... ! Package amsmath Error: \begin{split} won’t work here. See the amsmath package documentation for explanation. Type H <return> for immediate help. ... l.22 ...end{split}\notag\\\begin{split}\end{split} \notag ! Misplaced \cr. \math@cr@@@ ->\cr l.23 \end{gather} ! Misplaced \noalign. \math@cr@@ ... \iffalse }\fi \math@cr@@@ \noalign {\vskip #1\relax } l.23 \end{gather} ! Missing $ inserted. <inserted text> $ l.23 \end{gather} ! Missing } inserted. <inserted text> } l.23 \end{gather} ! Extra }, or forgotten \endgroup. \math@cr@@ ...th@cr@@@ \noalign {\vskip #1\relax } l.23 \end{gather} ! Misplaced \noalign. \black@ #1->\noalign {\ifdim #1>\displaywidth \dimen@ \prevdepth \nointerlin... l.23 \end{gather} ! Extra }, or forgotten \endgroup. \endgather ->\math@cr \black@ \totwidth@ \egroup $$\ignorespacesafterend l.23 \end{gather} ! LaTeX Error: \begin{preview} on input line 15 ended by \end{gather}. See the LaTeX manual or LaTeX Companion for explanation. Type H <return> for immediate help. ... l.23 \end{gather} ! Missing $ inserted. <inserted text> $ l.23 \end{gather} ! Display math should end with $$. <to be read again> \endgroup l.23 \end{gather} [1] ! LaTeX Error: \begin{document} ended by \end{preview}. See the LaTeX manual or LaTeX Companion for explanation. Type H <return> for immediate help. ... l.24 \end{preview} ! Extra \endgroup. <recently read> \endgroup l.24 \end{preview} (./math.aux) ) (see the transcript file for additional information) Output written on math.dvi (1 page, 144 bytes). Transcript written on math.log.

Examples

>>> V = PoschlTeller(v0=-1.0,r0=2.0/6.2527)
>>> V.gE(E=1.0)
0.1599309...
>>> V.k_cot_delta_analytic(E=1.0)
0.1599309...
mu0[source]

Stefano, Alex, etc. use this parameter instead, so we provide an alias.

re_ = <numpy.lib.function_base.vectorize object at 0x9975950>[source]
set_ainv(ainv=0.0)[source]

Set the coefficient v0 so that the inverse scattering length is ainv.

>>> V = PoschlTeller()
>>> V.v0 = -0.5
>>> V.set_ainv(0.0)
>>> np.allclose(V.v0, -1.0)
True
set_ainv_re(ainv=0.0, re=0.07)[source]

Set the coefficients v0 and r0 so that the inverse scattering length is ainv and the effective range is re.

>>> V = PoschlTeller(r0=1.0)
>>> re = 1.0
>>> V.set_ainv_re(ainv=0.0,re=re)
>>> np.allclose(V.r0, re)
True
>>> np.allclose(V.v0, -1.0)
True
class mmf.physics.potentials.two_body.BargmannPotential(a=1.0, re=0.4, hbar=0.15915494309189535, M=0.5, **kwarg)[source]

Bases: mmf.physics.potentials.two_body.Potential

Bargmann Potential where the effective range expansion is exact with only a and r_e.

>> V = BargmannPotential(a=1.0,re=0.01) >> (a,re) = V.get_ar() >> (a,re) (1.0, 0.1)

Methods

__call__(r)
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V’(r)/V(r).
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.
__init__(a=1.0, re=0.4, hbar=0.15915494309189535, M=0.5, **kwarg)[source]
dVdV(r)[source]

Return V’(r)/V(r).

V:=v0*beta*exp(-2*b*r)/(1+beta*exp(-2*b*r))^2; with(codegen); C(simplify(diff(V,r)/V),optimized);

class mmf.physics.potentials.two_body.PotentialV0R0(abs_tol=1e-12, rel_tol=1e-12, Rmin=None)[source]

Bases: mmf.physics.potentials.two_body.Potential

Mixing for two-parameter potentials with a parameter v_0 governing the overall strength and a parameter r_0 governing the effective range. The potential can have additional parameters, but these two will be adjusted to set the scattering-length and range.

Methods

__call__(r) Return V(r).
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V’(r)/V(r). Should be overloaded for speed and/or
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.
set_ainv_re(ainv, re) Uses a non-linear root-finder to set the parameters v0 and r0 so that the scattering length is 1.0/ainv and the effective range is re.

Construct the potential object.

Methods

__call__(r) Return V(r).
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V’(r)/V(r). Should be overloaded for speed and/or
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.
set_ainv_re(ainv, re) Uses a non-linear root-finder to set the parameters v0 and r0 so that the scattering length is 1.0/ainv and the effective range is re.
__init__(abs_tol=1e-12, rel_tol=1e-12, Rmin=None)

Construct the potential object.

set_ainv_re(ainv, re)[source]

Uses a non-linear root-finder to set the parameters v0 and r0 so that the scattering length is 1.0/ainv and the effective range is re.

>>> V = PoschlTeller()
>>> V.set_ainv_re(ainv=0.0,re=1.0)
>>> np.allclose(V.r0, 1.0)
True
>>> np.allclose(V.v0, -1.0)
True
class mmf.physics.potentials.two_body.DoubleGaussian(v0=0.7860976072874742, r0=0.6666666666666666, mu0=None, f_v=-4.0, f_r=2.0, hbar=1.0, M=0.5, abs_tol=1e-12, rel_tol=1e-12, Rmin=0.1)[source]

Bases: mmf.physics.potentials.two_body.PotentialV0R0

Double Gaussian potential:

V(r) &= -v_0 \frac{\hbar^2}{M} \mu_0^2 \left(
  \exp[-(\mu_0 r/2)^{2}] + f_v \exp[-(\mu_0 r f_r/2)^2]
  \right)\\
    &= -v_0 \frac{\hbar^2}{M} \mu_0^2 \left(
  \exp[-(\mu_0 r/2)^{2}] + f_v \exp[-(\mu_0 r f_r/2)^2]
  \right)\\

Note that \mu_0 = 4/r_0 here.

>>> v = DoubleGaussian(mu0=12.0)
>>> v.get_ar()              # One of Stefano's values to check.
(-2887043864..., 0.3293...)
>>> v.set_ainv_re(ainv=0.0, re=1.0)
>>> v.v0, v.r0
(0.786097607..., 1.012180386...)

Methods

__call__(r)
compute_R() Return R >= self.Rmin such that `abs(self(R)) <
dVdV(r) Return V’(r)/V(r). Should be overloaded for speed and/or
gE(E[, R, method]) Return k*cot(delta) as a function of energy E for the specified
get_ar([R, method]) Return (a,re) where a is the s-wave scattering length and re is the s-wave effective range.
get_deq(f, jac, U0[, method]) Return a “deq” object capable of integrating the system
get_deq_E(E[, method]) Return solution to energy dependent differential equation for (u,du) for finding the s-wave phase shift.
get_deq_ar([method]) Return the differential equation (u,du,uu,duu) for finding the s-wave scattering length and the s-wave effective range.
get_deq_gE(E[, method]) Return solution to energy dependent differential equation
get_deq_gE_old(E[, method, version]) Old version in terms of fp and fm.
get_mu0()
plot([Rmax]) Plot potential.
plot_gE(Erange[, Rmax]) Plot potential.
plot_u(E[, Rmax]) Plot radial wavefunction.
set_ainv_re(ainv, re) Uses a non-linear root-finder to set the parameters v0 and r0 so that the scattering length is 1.0/ainv and the effective range is re.
set_mu0(mu0)
__init__(v0=0.7860976072874742, r0=0.6666666666666666, mu0=None, f_v=-4.0, f_r=2.0, hbar=1.0, M=0.5, abs_tol=1e-12, rel_tol=1e-12, Rmin=0.1)[source]
get_mu0()[source]
mu0
set_mu0(mu0)[source]
mmf.physics.potentials.two_body.SoonYongPotential

alias of PoschlTeller