AtomicOrbital

class sisl.AtomicOrbital(*args, **kwargs)[source]

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.

Parameters
*argslist of arguments

list of arguments can be in different input options

q0float, optional

initial charge

tagstr, optional

user defined tag

Examples

>>> r = np.linspace(0, 5, 50)
>>> f = np.exp(-r)
>>> #                    n, l, m, [Z, [P]]
>>> orb1 = AtomicOrbital(2, 1, 0, 1, (r, f))
>>> orb2 = AtomicOrbital(n=2, l=1, m=0, Z=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
Attributes
Rfloat

maximum radius (in Ang)

nint

principal quantum number

lint

azimuthal quantum number

mint

magnetic quantum number

Zint

zeta shell

Pbool

whether this corresponds to a polarized shell (True)

ffunc

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

q0float

initial electronic charge

tagstr

user defined tag

Attributes

P

R

Z

f

l

m

n

orb

q0

tag

Methods

__init__(self, \*args, \*\*kwargs)

Initialize atomic orbital object

copy(self)

Create an exact copy of this object

equal(self, other[, psi, radial])

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

name(self[, tex])

Return named specification of the atomic orbital

psi(self, r)

Calculate \(\phi(\mathbf r)\) at a given point (or more points)

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

Calculate \(\phi(|\mathbf R|, \theta, \phi)\) at a given point (in spherical coordinates)

radial(self, r[, is_radius])

Calculate the radial part of the wavefunction \(f(\mathbf R)\)

scale(self, scale)

Scale the orbital by extending R by scale

set_radial(self, \*args)

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

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

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

toGrid(self[, precision, c, R, dtype, Z])

Create a Grid with only this orbital wavefunction on it

toSphere(self[, center])

Return a sphere with radius equal to the orbital size

P
R
Z
copy(self)[source]

Create an exact copy of this object

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

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

Parameters
otherOrbital

comparison orbital

psibool, optional

also compare that the full psi are the same

radialbool, optional

also compare that the radial parts are the same

property f
l
m
n
name(self, tex=False)[source]

Return named specification of the atomic orbital

orb
psi(self, r)[source]

Calculate \(\phi(\mathbf r)\) at a given point (or more points)

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

Parameters
rarray_like

the vector from the orbital origin

Returns
numpy.ndarray

basis function value at point r

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

Calculate \(\phi(|\mathbf R|, \theta, \phi)\) at a given point (in spherical coordinates)

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

Parameters
rarray_like

the radius from the orbital origin

thetaarray_like

azimuthal angle in the \(x-y\) plane (from \(x\))

phiarray_like

polar angle from \(z\) axis

cos_phibool, optional

whether phi is actually \(cos(\phi)\) which will be faster because cos is not necessary to call.

Returns
numpy.ndarray

basis function value at point r

q0
radial(self, r, is_radius=True)[source]

Calculate the radial part of the wavefunction \(f(\mathbf R)\)

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

Parameters
rarray_like

radius from the orbital origin, for is_radius=False r must be vectors

is_radiusbool, optional

whether r is a vector or the radius

Returns
numpy.ndarray

radial orbital value at point r

scale(self, scale)

Scale the orbital by extending R by scale

set_radial(self, *args)[source]

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

See SphericalOrbital.set_radial where these arguments are passed to.

spher(self, theta, phi, cos_phi=False)[source]

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

Parameters
thetaarray_like

azimuthal angle in the \(x-y\) plane (from \(x\))

phiarray_like

polar angle from \(z\) axis

cos_phibool, optional

whether phi is actually \(cos(\phi)\) which will be faster because cos is not necessary to call.

Returns
numpy.ndarray

spherical harmonics at angles \(\theta\) and \(\phi\)

tag
toGrid(self, precision=0.05, c=1.0, R=None, dtype=<class 'numpy.float64'>, Z=1)

Create a Grid with only this orbital wavefunction on it

Parameters
precisionfloat, optional

used separation in the Grid between voxels (in Ang)

cfloat or complex, optional

coefficient for the orbital

Rfloat, optional

box size of the grid (default to the orbital range)

dtypenumpy.dtype, optional

the used separation in the Grid between voxels

Zint, optional

atomic number associated with the grid

toSphere(self, center=None)

Return a sphere with radius equal to the orbital size

Returns
~sisl.shape.Sphere

sphere with a radius equal to the radius of this orbital