sisl.AtomicOrbital

class sisl.AtomicOrbital

Bases: Orbital

A projected atomic orbital consisting of real harmonics

The AtomicOrbital is a specification of the SphericalOrbital by assigning the magnetic quantum number m to the object.

AtomicOrbital should always be preferred over the SphericalOrbital because it explicitly contains all quantum numbers.

The atomic orbital has a radial part defined by an external function; this is then expanded using spherical harmonics

Ylm(θ,φ)=(1)m2l+14π(lm)!(l+m)!eimθPlm(cos(φ))ϕlmn(r)=R(|r|)Ylm(θ,φ)

where the function R(|r|) is user-defined.

Parameters:
  • *args (list of arguments) – list of arguments can be in different input options

  • R – See Orbital for details.

  • q0 – initial charge

  • tag – user defined tag

Examples

>>> r = np.linspace(0, 5, 50)
>>> f = np.exp(-r)
>>> #                    n, l, m, [zeta, [P]]
>>> orb1 = AtomicOrbital(2, 1, 0, 1, (r, f))
>>> orb2 = AtomicOrbital(n=2, l=1, m=0, zeta=1, (r, f))
>>> orb3 = AtomicOrbital("2pzZ", (r, f))
>>> orb4 = AtomicOrbital("2pzZ1", (r, f))
>>> orb5 = AtomicOrbital("pz", (r, f))
>>> orb2 == orb3
True
>>> orb2 == orb4
True
>>> orb2 == orb5
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 named specification of the atomic orbital

psi(r)

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

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

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

radial(r, *args, **kwargs)

Calculate the radial part of the wavefunction f(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[, cos_phi])

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

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

P

Whether this is polarized shell or not

R

Maxmimum radius of orbital

l

l quantum number

m

m quantum number

n

n shell

orb

Orbital with radial part

q0

Initial charge

tag

Named tag of orbital

zeta

ζ shell

copy()

Create an exact copy of this object

Parameters:

orbital (AtomicOrbital)

Return type:

AtomicOrbital

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) – also compare that the full psi are the same

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

name(tex=False)[source]

Return named specification of the atomic orbital

psi(r)[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) – the vector from the orbital origin

Returns:

basis function value at point r

Return type:

numpy.ndarray

psi_spher(r, theta, phi, 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

  • 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)[source]

Calculate the radial part of the wavefunction f(r)

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

Parameters:

r (ndarray) – radius from the orbital origin

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)[source]

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

See SphericalOrbital.set_radial where these arguments are passed to.

spher(theta, phi, 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

  • 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 ϕ

Return type:

numpy.ndarray

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 P

Whether this is polarized shell or not

property R

Maxmimum radius of orbital

property l

l quantum number

property m

m quantum number

property n

n shell

property orb

Orbital with radial part

property q0

Initial charge

property tag

Named tag of orbital

property zeta

ζ shell