SphericalOrbital¶
-
class
sisl.
SphericalOrbital
(l, rf_or_func, q0=0.0, tag='')[source]¶ An arbitrary orbital class where \(\phi(\mathbf r)=f(|\mathbf r|)Y_l^m(\theta,\varphi)\)
Note that in this case the used spherical harmonics is:
\[Y^m_l(\theta,\varphi) = (-1)^m\sqrt{\frac{2l+1}{4\pi} \frac{(l-m)!}{(l+m)!}} e^{i m \theta} P^m_l(\cos(\varphi))\]The resulting orbital is
\[\phi_{lmn}(\mathbf r) = f(|\mathbf r|) Y^m_l(\theta, \varphi)\]where typically \(f(|\mathbf r|)\equiv\phi_{ln}(|\mathbf 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 \(\phi_{ln}(|\mathbf r|)\) function.
- Parameters
- lint
azimuthal quantum number
- rf_or_functuple of (r, f) or func
radial components as a tuple/list, or the function which can interpolate to any R See
set_radial
for details.- q0float, optional
initial charge
- tagstr, optional
user defined tag
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
- Attributes
- Rfloat
maximum radius (in Ang)
- q0float
initial electronic charge
- lint
azimuthal quantum number
- ffunc
interpolation function that returns f(r) for a given radius
- tagstr
user defined tag
Attributes
Methods
__init__
(self, l, rf_or_func[, q0, tag])Initialize spherical 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 a named specification of the orbital (
tag
)psi
(self, r[, m])Calculate \(\phi(\mathbf R)\) at a given point (or more points)
psi_spher
(self, r, theta, phi[, m, 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, \*\*kwargs)Update the internal radial function used as a \(f(|\mathbf r|)\)
spher
(self, theta, phi[, m, cos_phi])Calculate the spherical harmonics of this orbital at a given point (in spherical coordinates)
toAtomicOrbital
(self[, m, n, Z, P, q0])Create a list of
AtomicOrbital
objectstoGrid
(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
-
R
¶
-
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
-
f
¶
-
l
¶
-
psi
(self, r, m=0)[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 of (:, 3)
the vector from the orbital origin
- mint, optional
magnetic quantum number, must be in range
-self.l <= m <= self.l
- Returns
- numpy.ndarray
basis function value at point r
-
psi_spher
(self, r, theta, phi, m=0, 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
- mint, optional
magnetic quantum number, must be in range
-self.l <= m <= self.l
- 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
-
set_radial
(self, *args, **kwargs)[source]¶ Update the internal radial function used as a \(f(|\mathbf 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 off(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, fnumpy.ndarray
the radial positions and the radial function values at r.
- funccallable
a function which enables evaluation of the radial function. The function should accept a single array and return a single array.
- interpcallable, 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.
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.f(R) >>> o.set_radial(r, f, interp=i_interp1d) >>> f_interp1d = o.f(R) >>> o.set_radial(r, f, interp=i_spline) >>> f_spline = o.f(R) >>> np.allclose(f_univariate, f_interp1d) True >>> np.allclose(f_univariate, f_spline) True
-
spher
(self, theta, phi, m=0, 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
- mint, optional
magnetic quantum number, must be in range
-self.l <= m <= self.l
- 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\) and given quantum number m
-
tag
¶
-
toAtomicOrbital
(self, m=None, n=None, Z=1, P=False, q0=None)[source]¶ Create a list of
AtomicOrbital
objectsThis 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
- mint or list or None
if
None
it defaults to-l:l
, else only for the requested m- Zint, optional
the specified zeta-shell
- Pbool, optional
whether the orbitals are polarized.
- q0float, optional
the initial charge per orbital, initially \(q_0 / (2l+1)\) with \(q_0\) from this object
- Returns
- AtomicOrbitalfor passed m an atomic orbital will be returned
- list of AtomicOrbitalfor each \(m\in[-l;l]\) an atomic orbital will be returned in the list
-
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