sisl.SphericalOrbital

class sisl.SphericalOrbital

Bases: Orbital

An arbitrary orbital class which only contains the harmonical part of the wavefunction where ϕ(r)=f(|r|)Ylm(θ,φ)

Note that in this case the used spherical harmonics is:

Ylm(θ,φ)=(1)m2l+14π(lm)!(l+m)!eimθPlm(cos(φ))

The resulting orbital is

ϕlmn(r)=f(|r|)Ylm(θ,φ)

where typically f(|r|)ϕln(|r|). The above equation clarifies that this class is only intended for each l, and that subsequent m orders may be extracted by altering the spherical harmonic. Also, the quantum number n is not necessary as that value is implicit in the ϕln(|r|) function.

Parameters:
  • l – azimuthal quantum number

  • rf_or_func – radial components as a tuple/list, or the function which can interpolate to any R See set_radial for details.

  • R – See Orbital for details.

  • q0 – initial charge

  • tag – user defined tag

  • **kwargs – arguments passed directly to set_radial(rf_or_func, **kwargs)

f

interpolation function that returns f(r) for a given radius

Type:

func

Examples

>>> from scipy.interpolate import interp1d
>>> orb = SphericalOrbital(1, (np.arange(10.), np.arange(10.)))
>>> orb.equal(SphericalOrbital(1, interp1d(np.arange(10.), np.arange(10.),
...       fill_value=(0., 0.), kind="cubic", bounds_error=False)))
True

Methods

copy()

Create an exact copy of this object

equal(other[, psi, radial])

Compare two orbitals by comparing their radius, and possibly the radial and psi functions

name([tex])

Return a named specification of the orbital (tag)

psi(r[, m])

Calculate ϕ(r) at a given point (or more points)

psi_spher(r, theta, phi[, m, cos_phi])

Calculate ϕ(|r|,θ,ϕ) at a given point (in spherical coordinates)

radial(r, *args, **kwargs)

Calculate the radial part of spherical orbital R(r)

scale(scale)

Scale the orbital by extending R by scale

set_radial(*args, **kwargs)

Update the internal radial function used as a f(|r|)

spher(theta, phi[, m, cos_phi])

Calculate the spherical harmonics of this orbital at a given point (in spherical coordinates)

toAtomicOrbital([m, n, zeta, P, q0])

Create a list of AtomicOrbital objects

toGrid([precision, c, R, dtype, atom])

Create a Grid with only this orbital wavefunction on it

toSphere([center])

Return a sphere with radius equal to the orbital size

Attributes

R

Maxmimum radius of orbital

l

l quantum number

q0

Initial charge

tag

Named tag of orbital

copy()

Create an exact copy of this object

Parameters:

orbital (SphericalOrbital)

Return type:

SphericalOrbital

equal(other, psi=False, radial=False)[source]

Compare two orbitals by comparing their radius, and possibly the radial and psi functions

Parameters:
  • other (Orbital) – comparison orbital

  • psi (bool, optional) – also compare that the full psi are the same

  • radial (bool, optional) – also compare that the radial parts are the same

name(tex=False)

Return a named specification of the orbital (tag)

psi(r, m=0)[source]

Calculate ϕ(r) at a given point (or more points)

The position r is a vector from the origin of this orbital.

Parameters:
  • r (ndarray of (:, 3)) – vector from the orbital origin

  • m (int) – magnetic quantum number, must be in range -self.l <= m <= self.l

Returns:

basis function value at point r

Return type:

numpy.ndarray

psi_spher(r, theta, phi, m=0, cos_phi=False)[source]

Calculate ϕ(|r|,θ,ϕ) at a given point (in spherical coordinates)

This is equivalent to psi however, the input is given in spherical coordinates.

Parameters:
  • r (ndarray) – the radius from the orbital origin

  • theta (ndarray) – azimuthal angle in the xy plane (from x)

  • phi (ndarray) – polar angle from z axis

  • m (int) – magnetic quantum number, must be in range -self.l <= m <= self.l

  • cos_phi (bool) – whether phi is actually cos(ϕ) which will be faster because cos is not necessary to call.

Returns:

basis function value at point r

Return type:

numpy.ndarray

radial(r, *args, **kwargs)

Calculate the radial part of spherical orbital R(r)

The position r is a vector from the origin of this orbital.

Parameters:
  • r (ndarray) – radius from the orbital origin

  • *args – arguments passed to the radial function

  • **args – keyword arguments passed to the radial function

Returns:

radial orbital value at point r

Return type:

numpy.ndarray

scale(scale)

Scale the orbital by extending R by scale

Parameters:
Return type:

Orbital

set_radial(*args, **kwargs)

Update the internal radial function used as a f(|r|)

This can be called in several ways:

set_radial(r, f)

which uses scipy.interpolate.UnivariateSpline(r, f, k=3, s=0, ext=1, check_finite=False) to define the interpolation function (see interp keyword). Here the maximum radius of the orbital is the maximum r value, regardless of f(r) is zero for smaller r.

set_radial(func)

which sets the interpolation function directly. The maximum orbital range is determined automatically to a precision of 0.0001 AA.

Parameters:
  • r (numpy.ndarray) – the radial positions and the radial function values at r.

  • f (numpy.ndarray) – the radial positions and the radial function values at r.

  • func (callable) – a function which enables evaluation of the radial function. The function should accept a single array and return a single array.

  • interp (callable, optional) – When two non-keyword arguments are passed this keyword will be used. It is the interpolation function which should return the equivalent of func. By using this one can define a custom interpolation routine. It should accept two arguments, interp(r, f) and return a callable that returns interpolation values. See examples for different interpolation routines.

Return type:

None

Examples

>>> from scipy import interpolate as interp
>>> o = SphericalOrbital(1, lambda x:x)
>>> r = np.linspace(0, 4, 300)
>>> f = np.exp(-r)
>>> def i_univariate(r, f):
    ...    return interp.UnivariateSpline(r, f, k=3, s=0, ext=1, check_finite=False)
>>> def i_interp1d(r, f):
    ...    return interp.interp1d(r, f, kind="cubic", fill_value=(f[0], 0.), bounds_error=False)
>>> def i_spline(r, f):
    ...    from functools import partial
...    tck = interp.splrep(r, f, k=3, s=0)
...    return partial(interp.splev, tck=tck, der=0, ext=1)
>>> R = np.linspace(0, 4, 400)
>>> o.set_radial(r, f, interp=i_univariate)
>>> f_univariate = o.radial(R)
>>> o.set_radial(r, f, interp=i_interp1d)
>>> f_interp1d = o.radial(R)
>>> o.set_radial(r, f, interp=i_spline)
>>> f_spline = o.radial(R)
>>> np.allclose(f_univariate, f_interp1d)
True
>>> np.allclose(f_univariate, f_spline)
True
spher(theta, phi, m=0, cos_phi=False)[source]

Calculate the spherical harmonics of this orbital at a given point (in spherical coordinates)

Parameters:
  • theta (ndarray) – azimuthal angle in the xy plane (from x)

  • phi (ndarray) – polar angle from z axis

  • m (int) – magnetic quantum number, must be in range -self.l <= m <= self.l

  • cos_phi (bool) – whether phi is actually cos(ϕ) which will be faster because cos is not necessary to call.

Returns:

spherical harmonics at angles θ and ϕ and given quantum number m

Return type:

numpy.ndarray

toAtomicOrbital(m=None, n=None, zeta=1, P=False, q0=None)[source]

Create a list of AtomicOrbital objects

This defaults to create a list of AtomicOrbital objects for every m (for m in -l:l). One may optionally specify the sub-set of m to retrieve.

Parameters:
  • m (int or list or None) – if None it defaults to -l:l, else only for the requested m

  • zeta (int) – the specified ζ-shell

  • n (int | None) – specify the n quantum number

  • P (bool) – whether the orbitals are polarized.

  • q0 (float | None) – the initial charge per orbital, initially q0/(2l+1) with q0 from this object

Returns:

  • AtomicOrbital (for passed `m an atomic orbital will be returned`)

  • list of AtomicOrbital (for each :math:`min[-l;l] an atomic orbital will be returned in the list`)

toGrid(precision=0.05, c=1.0, R=None, dtype=np.float64, atom=1)

Create a Grid with only this orbital wavefunction on it

Parameters:
  • precision (float, optional) – used separation in the Grid between voxels (in Ang)

  • c (float or complex, optional) – coefficient for the orbital

  • R (float, optional) – box size of the grid (default to the orbital range)

  • dtype (numpy.dtype, optional) – the used separation in the Grid between voxels

  • atom (optional) – atom associated with the grid; either an atom instance or something that Atom(atom) would convert to a proper atom.

toSphere(center=None)

Return a sphere with radius equal to the orbital size

Returns:

sphere with a radius equal to the radius of this orbital

Return type:

Sphere

property R

Maxmimum radius of orbital

property l

l quantum number

property q0

Initial charge

property tag

Named tag of orbital