moldesign.orbitals package¶
-
class
moldesign.orbitals.
BasisSet
(mol, orbitals, name=None, h1e=None, overlaps=None, angulartype=None, **kwargs)[source]¶ Bases:
moldesign.orbitals.orbitals.MolecularOrbitals
Stores a basis, typically of atomic orbitals.
This is a special orbital type
-
density_matrix
¶
-
fock
¶
-
Submodules¶
moldesign.orbitals.basis module¶
-
class
moldesign.orbitals.basis.
BasisSet
(mol, orbitals, name=None, h1e=None, overlaps=None, angulartype=None, **kwargs)[source]¶ Bases:
moldesign.orbitals.orbitals.MolecularOrbitals
Stores a basis, typically of atomic orbitals.
This is a special orbital type
-
density_matrix
¶
-
fock
¶
-
moldesign.orbitals.gaussians module¶
Note: this code is currently unused and untested and will be refactored soon
-
class
moldesign.orbitals.gaussians.
AbstractFunction
[source]¶ Bases:
object
Abstract base class for basis functions
-
overlap
(other, normalized=False)[source]¶ Overlap of this function with another:
\[\int f_1(\mathbf r) f_2(\mathbf r) d^N \mathbf r\]Parameters: - other (AbstractFunction) –
- normalized (bool) – If True, return the overlap of the two NORMALIZED functions.
Returns: value of the overlap
Return type:
-
-
class
moldesign.orbitals.gaussians.
AtomicBasisFunction
(atom, n=None, l=None, m=None, cart=None, primitives=None)[source]¶ Bases:
moldesign.orbitals.orbitals.Orbital
-
aotype
¶ A string describing the orbital’s state.
Examples
>>> AtomicBasisFunction(n=1, l=0).aotype '1s' >>> AtomicBasisFunction(n=2, l=1, cart='y').aotype '2py' >>> AtomicBasisFunction(n=3, l=2, m=0).aotype '3d(z^2)'
-
name
¶
-
num_primitives
¶
-
orbtype
¶ A string describing the orbital’s angular momentum state.
Examples
>>> AtomicBasisFunction(n=1, l=0).orbtype 's' >>> AtomicBasisFunction(n=2, l=1, cart='y').orbtype 'py' >>> AtomicBasisFunction(n=3, l=2, m=0).orbtype 'd(z^2)'
-
-
class
moldesign.orbitals.gaussians.
CartesianGaussian
(center, exp, powers, coeff=None)[source]¶ Bases:
moldesign.orbitals.gaussians.AbstractFunction
Stores an N-dimensional gaussian function of the form:
\[G(\mathbf r) = C \times \left( \prod_{i=1}^N{{r_i}^{p_i} } \right) e^{-a |\mathbf r - \mathbf{r}_0|^2}\]For a three-dimensional gaussian, this is
- ..math::
- G(x,y,z) = C times x^{p_1} y^{p_2} z^{p_3} e^{-a |mathbf r - mathbf{r}_0|^2}
where C is
self.coeff
, a isself.exp
, r0 isself.center
, and \(p_1, p_2, ...\) are given in the arrayself.powers
References
Levine, Ira N. Quantum Chemistry, 5th ed. Prentice Hall, 2000. 486-94.
-
angular
¶ Angular momentum of this function (sum of cartesian powers)
-
integral
¶ Integral of this this gaussian over all N-dimensional space.
This is implemented only for 0 and positive integer cartesian powers. The integral is 0 if any of the powers are odd. Otherwise, the integral is given by: .. math:
\int G(\mathbf r) d^N \mathbf r & = c \int d^N e^{-a x^2} \mathbf r \prod_{i=1}^N{{r_i}^{p_i} } \\ &= (2a)^{-\sum_i p_i} \left( \frac{\pi}{2 a} \right) ^ {N/2} \prod_{i=1}^N{(p_i-1)!!}
where N is the dimensionality of the gaussian, \(p_i\) are the cartesian powers, and _!!_ is the “odd factorial” (\(n!!=1\times 3\times 5 \times ... \times n\))
References
- Dwight, Herbert B. Tables of Integrals and other Mathematical Data, 3rd ed.
- Macmillan 1957. 201.
-
ndim
¶
-
ndims
¶
-
num_dimensions
¶
-
moldesign.orbitals.gaussians.
Gaussian
(center, exp, coeff=1.0)[source]¶ Constructor for an N-dimensional gaussian function.
The function is given by: .. math:
G(\mathbf r) = C e^{-a\left| \mathbf r - \mathbf r_0 \right|^2}
where C is
self.coeff
, a isself.exp
, and \(\mathbf r_0\) isself.center
.Note
This is just a special case of a cartesian gaussian where all the powers are 0.
-
class
moldesign.orbitals.gaussians.
SphericalGaussian
(center, exp, n, l, m, coeff=None)[source]¶ Bases:
moldesign.orbitals.gaussians.AbstractFunction
Stores a 3-dimensional spherical gaussian function:
\[G_{nlm}(\mathbf r) = C Y^l_m(\mathbf r - \mathbf r_0) r^n e^{-a\left| \mathbf r - \mathbf r_0 \right|^2}\]where C is
self.coeff
, a isself.exp
, and \(\mathbf r_0\) isself.center
, (n,l,m) are given by(self.n, self.l, self.m)
, and \(Y^l_m(\mathbf{r})\) is a spherical harmonic.Note
self.SPHERE_TO_CART stores expansion coefficients for spherical gaussians in terms of cartesian gaussians. They are taken from the reference below.
References
- Schlegel and Frisch. Transformation between cartesian and pure spherical harmonic gaussians.
- Int J Quantum Chem 54, 83-87 (1995). doi:10.1002/qua.560540202
moldesign.orbitals.orbitals module¶
Class definitions for atomic and molecular orbitals.
Notes
- In this documentation, we use the following conventions for labeling orbitals:
- atomic orbitals using lower case greek labels and subscripts, e.g.,
:math:`left| mu
ight angle, F_{ u lambda}, etc.
- molecular orbitals use lower case labels and subscripts, e.g.,
:math:`left| i
ight angle, F_{kl}, etc.
- adiabatic electronic states are indexed using capital letters, _N_, _L_, _M_, etc.
-
class
moldesign.orbitals.orbitals.
MolecularOrbitals
(orbitals, wfn=None, basis=None, canonical=False, orbtype=None)[source]¶ Bases:
object
Stores a wfn of molecular orbitals in an AO wfn Orbitals are accessed as orbs[orbital index, ao index]
-
align_phases
(other, threshold=0.5, assert_same_type=True)[source]¶ Flip the signs of these orbitals to bring them into maximum coincidence with another set of orbitals
Parameters: - other (MolecularOrbitals) – the “reference” set of orbitals to match phases with
- threshold (float) – only flip orbital if the overlap is less than -1*threshold
- assert_same_type (bool) – require that
self.orbtype == other.orbtype
Note
This function assumes that the overlap matrix is the same for both sets of orbitals - this is a reasonable assumption if the two sets of orbitals were calculated at very similar molecular geometries.
-
energies
¶ u.Vector[energy] – energies of the molecular orbitals
This is just the diagonal of the fock matrix
-
fock
¶ u.Array[energy] – Fock matrix for these orbitals
-
from_ao
(ao_operator)[source]¶ Transform an operator into this orbital basis from the ao basis
Given the matrix elements :math:`hat O_{mu- u}` of an operator over AO basis indices
- :math:`mu,
- u`, returns the operator’s matrix elements \(\hat O_{ij}\) over
orbital indices \(i,j\):
- ..math::
- hat O_{ij} =
- left langle i
ight| hat O left| j ight angle =
sum_{muu}C_{i mu} O_{mu u} C_{j u}
where \(C_{i \mu}\) is the expansion coefficient for AO basis function \(\mu\) in molecular orbital _i_.
- Args:
- ao_operator (u.Array): matrix elements of the operator in the ao basis
- Returns:
- u.Array: matrix elements of the operator in this orbital basis
- Note:
- Assumes that this set of orbitals is orthogonal
-
h1e
¶ u.Array[energy] – 1-electron matrix elements for these orbitals
-
h2e
¶ u.Array[energy] – 2-electron matrix elements for these orbitals
-
occupations
¶ np.ndarray – orbital occupation numbers
-
overlap
(other)[source]¶ Calculate overlaps between this and another set of orbitals
Parameters: other (MolecularOrbitals) – Returns: overlaps between the two sets of orbitals Return type: numpy.ndarray Example
>>> canonical = mol.wfn.canonical >>> atomic = mol.wfn.basis >>> overlaps = canonical.overlap(atomic) >>> overlaps[i, j] == canonical.orbitals[i].overlap(atomic.orbitals[j]) True
-
overlaps
¶ np.array – overlap matrix for these orbitals
-
to_ao
(mo_operator)[source]¶ Transform an operator from this orbital basis into the AO basis
Given the matrix elements \(\hat O_{ij}\) of an operator over orbital basis indices \(i,j\), returns the operator’s matrix elements :math:`hat O_{mu- u}` over
- orbital indices :math:`mu,
u`:
- ..math::
- hat O_{mu
- u} =
- left langle mu
ight| hat O left| u ight angle =
sum_{i,j,lambda,kappa}S_{mu lambda} C_{i lambda} O_{ij} C_{j kappa} S_{kappau}
where :math:`S_{muu} = left langle mu | u ight angle` is the AO overlap matrix
and \(C_{i \mu}\) is the expansion coefficient for AO basis function \(\mu\) in molecular orbital _i_.
- Args:
- mo_operator (u.Array): matrix elements of the operator in this orbital basis
- Returns:
- u.Array: matrix elements of the operator in the AO basis
-
-
class
moldesign.orbitals.orbitals.
Orbital
(coeffs, basis=None, wfn=None, occupation=None, name='unnamed')[source]¶ Bases:
object
Stores a single orbital and its meta-data Generally wants to be part of a set of MolecularOrbitals
The orbital is defined as .. math:
\left| i \right \rangle = \sum_\mu c_{i \mu} \left| \mu \right \rangle
where the coefficients \(c_{i \mu}\) are stored in
self.coeffs
and the basis orbitals \(\left| \mu \right \rangle\) are stored atself.basis
-
energy
¶ u.Scalar[energy] – This orbital’s energy
Note
This is equivalent to self.fock(self)
-
moldesign.orbitals.wfn module¶
-
class
moldesign.orbitals.wfn.
ElectronicWfn
(mol, num_electrons, model=None, aobasis=None, fock_ao=None, positions=None, civectors=None, description=None, density_matrix_ao=None)[source]¶ Bases:
object
Stores the results of a quantum chemistry calculation.
This is necessarily pretty flexible, but generally stores an LCAO wfn and one or more sets of orbitals. Can also store CI vectors, etc.
These objects will usually be created by quantum chemical energy models.
Parameters: - mol (moldesign.Molecule) – Molecule this wavefunction belongs to
- num_electrons (int) – number of electrons in this wavefunction
- model (moldesign.models.base.EnergyModelBase) – The model this wavefunction was created with
- aobasis (moldesign.orbitals.BasisSet) – The basis functions for the enclosed orbitals
- nbasis (int) – number of AO basis functions
- fock_ao (moldesign.units.Array[energy]) – fock matrix in the AO basis
- positions (moldesign.units.Array[length]) – positions of the nuclei for this wfn
- civectors (np.ndarray) – CI vectors (if applicable)
- description (str) – text describing the wfn (e.g. ‘RHF/STO-3G’, ‘CAS(2,2)/SA3/6-31G**’)
- density_matrix_ao (np.ndarray) – density matrix in the ao basis
-
align_orbital_phases
(other, assert_same=True)[source]¶ Align this wavefunction’s orbitals to have the same phase as those in other. :type other: ElectronicWfn :param assert_same: raise an exception if the two wavefunctions do not have the same kinds of orbitals
-
molecular_orbitals
¶ A synonym for self.orbitals[‘canonical’], since this is usually what’s wanted