3. API
3.1. pyscf.ita.promolecule module
- class pyscf.ita.promolecule.ProMolecule(method, element_methods=None, pro_charge=None, pro_multiplicity=None)[source]
Bases:
object
Promolecule class.
A promolecule is an approximation of a molecule constructed from a linear combination of atomic and/or ionic species. Properties of this promolecule can be computed from those of the atomic and/or ionic species, depending on whether the property is an extensive one or an intensive one.
Examples:
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> mf = scf.HF(mol) >>> mf.kernel() >>> pro_charge, pro_multiplicity = {'H':0, 'O':0}, {'H':2, 'O':3} >>> promol = ProMolecule.build(mf, pro_charge, pro_multiplicity)
- build(method=None, pro_charge=None, pro_multiplicity=None)[source]
Initialize a instance.
- Parameters:
method (PyscfMethod) – Pyscf method.
pro_charge (Dict{str:int}, optional) – Dictionary of charge for each element, by default None.
pro_multiplicity (Dict{str:int}, optional) – Dictionary of multiplicity for each element, by default None.
- Returns:
Instance of ProMolecule class.
- Return type:
- electron_density(grids_coords, spin='ab', deriv=2)[source]
Compute the electron density of the promolecule at the desired points.
- Parameters:
grid_coords (np.ndarray((N, 3), dtype=float)) – Points at which to compute the density.
spin (('ab' | 'a' | 'b' | 'm')) – Type of density to compute; either total, alpha-spin, beta-spin, or magnetization density, by default ‘ab’.
deriv (int) – List of molecule and promolecule derivative+1 level, by default 2.
- Returns:
Instance of PartitionDensity class.
- Return type:
- get_atom(element)[source]
Get the Pyscf Mole object for given element.
- Parameters:
element (str) – Symbol of element.
- Returns:
Pyscf Mole object for free atom.
- Return type:
Mole
- get_method(element)[source]
Get the Pyscf method for given element.
- Parameters:
element (str) – Symbol of element.
- Returns:
Pyscf method.
- Return type:
PyscfMethod
- get_scf(element)[source]
Get the mean field Pyscf method for given element.
- Parameters:
element (str) – Symbol of element.
- Returns:
Mean field pyscf method.
- Return type:
PyscfMethod
- get_xc(element)[source]
Get the exchange correlation functional for given element.
- Parameters:
element (str) – Symbol of element.
- Returns:
Name of exchange correlation functional.
- Return type:
str
3.2. pyscf.ita.aim module
- class pyscf.ita.aim.Becke(mol=None, grids=None)[source]
Bases:
AIM
Becke’s fuzzy atom approach.
- sharing_function(mol=None, grids=None)[source]
Becke sharing funcition \(\omega_A(\mathbf{r})\) defined as:
\[\omega_A(\mathbf{r}) = \frac{V_A (\mathbf{r}) }{\sum_A V_A (\mathbf{r})}\]- Parameters:
mol (Mole, optional) – Pyscf Mole instance, by default None.
grids (Grid, optional) – Pyscf Grids instance, by default None.
- Returns:
omega – Sharing functions for a list of atoms.
- Return type:
List[np.ndarray((N,), dtype=float)]
- class pyscf.ita.aim.Hirshfeld(promoldens=None)[source]
Bases:
AIM
Hirshfeld’s stockholder approach.
Examples:
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> mf = scf.HF(mol) >>> mf.kernel() >>> grids = dft.Grids(mol) >>> grids.build() >>> pro_charge, pro_multiplicity = {'H':0, 'O':0}, {'H':2, 'O':3} >>> aim = Hirshfeld.build(mf, grids.coords, pro_charge, pro_multiplicity) >>> omega = aim.sharing_function()
- classmethod build(method, grids_coords, pro_charge=None, pro_multiplicity=None, spin='ab', deriv=0)[source]
Class method to build the class.
- Parameters:
method (PyscfMethod) – Pyscf scf method or post-scf method instance.
grids_coords (np.ndarray((N,), dtype=float)) – Grids coords of N points.
pro_charge (Dict{str:int}, optional) – Dictionary of charge for each element, by default None.
pro_multiplicity (Dict{str:int}, optional) – Dictionary of multiplicity for each element, by default None.
spin (('ab' | 'a' | 'b' | 'm')) – Type of density to compute; either total, alpha-spin, beta-spin, or magnetization density, by default ‘ab’.
deriv (int) – List of molecule and promolecule derivative+1 level, by default 0.
- Returns:
Instance of Hirshfeld class.
- Return type:
obj
- sharing_function(promoldens=None)[source]
Hirshfeld sharing funcition \(\omega_A(\mathbf{r})\) defined as:
\[\omega_A(\mathbf{r}) = \frac{\rho^0_A (\mathbf{r}) }{\sum_A \rho^0_A (\mathbf{r})}\]- Parameters:
promoldens (PartitionDensity) – PartitionDensity instance for promolecule, by default None.
- Returns:
omega – Sharing functions for a list of atoms.
- Return type:
List[np.ndarray((N,), dtype=float)]
3.3. pyscf.ita.dens module
- class pyscf.ita.dens.ElectronDensity(rhos=None)[source]
Bases:
ndarray
Class to generate the electron density of the molecule at the desired points.
Examples:
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> grids = dft.gen_grid.Grids(mol) >>> grids.build() >>> mf = scf.HF(mol) >>> mf.kernel() >>> ed = ita.dens.ElectronDenstiy.build(mf, grids)
- classmethod build(method, grids_coords, spin='ab', deriv=2)[source]
Compute the electron density of the molecule at the desired points.
- Parameters:
method (PyscfMethod) – Pyscf scf method or post-scf method instance.
grid_coords (np.ndarray((N, 3), dtype=float)) – Points at which to compute the density.
spin (('ab' | 'a' | 'b' | 'm')) – Type of density to compute; either total, alpha-spin, beta-spin, or magnetization density, by default ‘ab’.
deriv (int) – List of molecule and promolecule derivative+1 level, by default 2.
- Returns:
obj – Instance of ElectronDensity class.
- Return type:
- density(mask=False, threshold=1e-30)[source]
Electron density \(\rho\left(\mathbf{r}\right)\).
- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
rho – Electron density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- gradient(mask=False, threshold=1e-30)[source]
Gradient of electron density \(\nabla \rho\left(\mathbf{r}\right)\).
This is the first-order partial derivatives of electron density w.r.t. coordinate \(\mathbf{r} = \left(x\mathbf{i}, y\mathbf{j}, z\mathbf{k}\right)\),
\[\nabla\rho\left(\mathbf{r}\right) = \left(\frac{\partial}{\partial x}\mathbf{i}, \frac{\partial}{\partial y}\mathbf{j}, \frac{\partial}{\partial z}\mathbf{k}\right) \rho\left(\mathbf{r}\right)\]- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
rho_grad
- Return type:
np.ndarray((3,N), dtype=float)
- gradient_norm(mask=False, threshold=1e-30)[source]
Norm of the gradient of electron density.
\[\lvert \nabla \rho\left(\mathbf{r}\right) \rvert = \sqrt{ \left(\frac{\partial\rho\left(\mathbf{r}\right)}{\partial x}\right)^2 + \left(\frac{\partial\rho\left(\mathbf{r}\right)}{\partial y}\right)^2 + \left(\frac{\partial\rho\left(\mathbf{r}\right)}{\partial z}\right)^2 }\]- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
grad_rho_norm
- Return type:
np.ndarray((N,), dtype=float)
- laplacian(mask=False, threshold=1e-30)[source]
Laplacian of electron density \(\nabla ^2 \rho\left(\mathbf{r}\right)\).
This is defined as the trace of Hessian matrix of electron density which is equal to the sum of its \(\left(\lambda_1, \lambda_2, \lambda_3\right)\) eigen-values:
\[\nabla^2 \rho\left(\mathbf{r}\right) = \nabla\cdot\nabla\rho\left(\mathbf{r}\right) = \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial x^2} + \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial y^2} + \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial z^2} = \lambda_1 + \lambda_2 + \lambda_3\]- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
rho_lapl
- Return type:
np.ndarray((N,), dtype=float)
- static spin_reduction(tensor, spin)[source]
Unified the shape of mo_coeff to square matirx, mo_occ to vector according to spin type.
- Parameters:
tensor (np.ndarray) – Tensor of mo_coeff or mo_occ.
spin (('ab' | 'a' | 'b' | 'm'), default='ab') – Type of density to compute; either total, alpha-spin, beta-spin, or magnetization density.
- Returns:
tensor – Square matrix mo_coeff or vector mo_occ.
- Return type:
np.ndarray
- tau(mask=False, threshold=1e-30)[source]
Laplacian of electron density \(\nabla ^2 \rho\left(\mathbf{r}\right)\).
This is defined as the trace of Hessian matrix of electron density which is equal to the sum of its \(\left(\lambda_1, \lambda_2, \lambda_3\right)\) eigen-values:
\[\nabla^2 \rho\left(\mathbf{r}\right) = \nabla\cdot\nabla\rho\left(\mathbf{r}\right) = \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial x^2} + \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial y^2} + \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial z^2} = \lambda_1 + \lambda_2 + \lambda_3\]- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
rho_tau
- Return type:
np.ndarray((N,), dtype=float)
- class pyscf.ita.dens.PartitionDensity(rhos=None)[source]
Bases:
list
Class to generate the partition electron density of the molecule at the desired points.
Examples:
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> grids = dft.gen_grid.Grids(mol) >>> grids.build() >>> mf = scf.HF(mol) >>> mf.kernel() >>> pd = ita.dens.PartitionDenstiy.orbital(mf, grids)
- density(mask=False, threshold=1e-30)[source]
Electron density \(\rho\left(\mathbf{r}\right)\).
- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
rho
- Return type:
np.ndarray((N,), dtype=float)
- gradient(mask=False, threshold=1e-30)[source]
Gradient of electron density \(\nabla \rho\left(\mathbf{r}\right)\).
This is the first-order partial derivatives of electron density w.r.t. coordinate \(\mathbf{r} = \left(x\mathbf{i}, y\mathbf{j}, z\mathbf{k}\right)\),
\[\nabla\rho\left(\mathbf{r}\right) = \left(\frac{\partial}{\partial x}\mathbf{i}, \frac{\partial}{\partial y}\mathbf{j}, \frac{\partial}{\partial z}\mathbf{k}\right) \rho\left(\mathbf{r}\right)\]- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
rho_grad
- Return type:
np.ndarray((3,N), dtype=float)
- gradient_norm(mask=False, threshold=1e-30)[source]
Norm of the gradient of electron density.
\[\lvert \nabla \rho\left(\mathbf{r}\right) \rvert = \sqrt{ \left(\frac{\partial\rho\left(\mathbf{r}\right)}{\partial x}\right)^2 + \left(\frac{\partial\rho\left(\mathbf{r}\right)}{\partial y}\right)^2 + \left(\frac{\partial\rho\left(\mathbf{r}\right)}{\partial z}\right)^2 }\]- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
grad_norm
- Return type:
np.ndarray((N,), dtype=float)
- laplacian(mask=False, threshold=1e-30)[source]
Laplacian of electron density \(\nabla ^2 \rho\left(\mathbf{r}\right)\).
This is defined as the trace of Hessian matrix of electron density which is equal to the sum of its \(\left(\lambda_1, \lambda_2, \lambda_3\right)\) eigen-values:
\[\nabla^2 \rho\left(\mathbf{r}\right) = \nabla\cdot\nabla\rho\left(\mathbf{r}\right) = \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial x^2} + \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial y^2} + \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial z^2} = \lambda_1 + \lambda_2 + \lambda_3\]- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
rho_lapl
- Return type:
np.ndarray((N,), dtype=float)
- classmethod orbital(method, grids_coords, spin='ab', deriv=1)[source]
Compute the orbital electron density of the molecule at the desired points.
- Parameters:
method (PyscfMethod) – Pyscf scf method or post-scf method instance.
grid_coords (np.ndarray((N, 3), dtype=float)) – Points at which to compute the density.
spin (('ab' | 'a' | 'b' | 'm')) – Type of density to compute; either total, alpha-spin, beta-spin, or magnetization density, by default ‘ab’.
deriv (int) – List of molecule and promolecule derivative+1 level, by default 1.
- Returns:
obj – Instance of PartitionDensity class.
- Return type:
- tau(mask=False, threshold=1e-30)[source]
Laplacian of electron density \(\nabla ^2 \rho\left(\mathbf{r}\right)\).
This is defined as the trace of Hessian matrix of electron density which is equal to the sum of its \(\left(\lambda_1, \lambda_2, \lambda_3\right)\) eigen-values:
\[\nabla^2 \rho\left(\mathbf{r}\right) = \nabla\cdot\nabla\rho\left(\mathbf{r}\right) = \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial x^2} + \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial y^2} + \frac{\partial^2\rho\left(\mathbf{r}\right)}{\partial z^2} = \lambda_1 + \lambda_2 + \lambda_3\]- Parameters:
mask (Bool) – If mask the corresponding element of the associated array which is less than given the threshold.
threshold (float) – Threshold for array element to mask or fill.
- Returns:
rho_tau
- Return type:
np.ndarray((N,), dtype=float)
3.4. pyscf.ita.ked module
- class pyscf.ita.ked.KineticEnergyDensity(dens=None, orbdens=None)[source]
Bases:
object
The kinetic-energy density class.
- dirac(rho=None, omega=None)[source]
Dirac kinetic energy density defined as:
\[\tau_\text{D} = \tfrac{3}{10} \left(6 \pi^2 \right)^{2/3} \left(\frac{\rho\left(\mathbf{r}\right)}{2}\right)^{5/3}\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ked – Kinetic energy density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- general(a=0.0, rho=None, rho_lapl=None, omega=None)[source]
Compute general(ish) kinetic energy density defined as:
\[\tau_\text{G} (\mathbf{r}, \alpha) = \tau_\text{PD} (\mathbf{r}) + \frac{1}{4} (a - 1) \nabla^2 \rho (\mathbf{r})\]- Parameters:
a (float) – Value of parameter \(a\).
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
rho_lapl (np.ndarray((N,), dtype=float), optional) – Electron density laplacian on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ked – Kinetic energy density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- gradient_expansion(a=None, b=None, rho=None, rho_grad_norm=None, rho_lapl=None, omega=None)[source]
Gradient expansion approximation of kinetic energy density defined as:
\[\tau_\text{GEA} \left(\mathbf{r}\right) = \tau_\text{TF} \left(\mathbf{r}\right) + a* \tau_\text{W} \left(\mathbf{r}\right) + b* \nabla^2 \rho\left(\mathbf{r}\right)\]There is a special case of
gradient_expansion()
with \(a=\tfrac{1}{9}\) and \(b=\tfrac{1}{6}\).- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
rho_grad_norm (np.ndarray((N,), dtype=float), optional) – Electron density graident norm on grid of N points, by default None.
rho_lapl (np.ndarray((N,), dtype=float), optional) – Electron density laplacian on grid of N points, by default None.
a (float, optional) – Value of parameter \(a\), by default None.
b (float, optional) – Value of parameter \(b\), by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ked – Kinetic energy density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- positive_definite(tau=None, omega=None)[source]
Positive definite or Lagrangian kinetic energy density defined as:
\[\tau_\text{PD} (\mathbf{r}) = \frac{1}{2} \sum_i^N n_i \rvert \nabla \phi_i (\mathbf{r}) \lvert^2\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ked – Kinetic energy density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- single_particle(rho_lapl=None, orbrho=None, orbrho_grad_norm=None, omega=None)[source]
Single-particle kinetic energy density defined as:
\[\tau_\text{s} = \sum_i \frac{1}{8} \frac{\nabla \rho_i \nabla \rho_i}{\rho_i} -\frac{1}{8} \nabla^2 \rho\]where \(\rho_i\) are the Kohn-Sham orbital densities.
- Parameters:
rho_lapl (np.ndarray((N,), dtype=float), optional) – Electron density laplacian on grid of N points, by default None.
orbrho (np.ndarray((N,), dtype=float), optional) – Electron density of orbitals on grid of N points, by default None.
orbrho_grad_norm (np.ndarray((N,), dtype=float), optional) – Electron density graident norm of orbitals on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ked – Kinetic energy density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- thomas_fermi(rho=None, omega=None)[source]
Thomas-Fermi kinetic energy density defined as:
\[\tau_\text{TF} \left(\mathbf{r}\right) = \tfrac{3}{10} \left(6 \pi^2 \right)^{2/3} \left(\frac{\rho\left(\mathbf{r}\right)}{2}\right)^{5/3}\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ked – Kinetic energy density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- weizsacker(rho=None, rho_grad_norm=None, omega=None)[source]
Weizsacker kinetic energy density defined as:
\[\tau_\text{W} (\mathbf{r}) = \frac{1}{8} \frac{\lvert \nabla\rho (\mathbf{r}) \rvert^2}{\rho (\mathbf{r})}\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
rho_grad_norm (np.ndarray((N,), dtype=float), optional) – Electron density graident norm on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ked – Kinetic energy density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
3.5. pyscf.ita.itad module
- class pyscf.ita.itad.ItaDensity(dens=None, prodens=None, keds=None)[source]
Bases:
object
Information-Theoretic Approch (ITA) Density class.
- G1(rho=None, prorho=None, rho_lapl=None, omega=None)[source]
G1 density defined as:
\[g_1(\mathbf{r}) \equiv {}^r_F i^{\prime}(\mathbf{r}) = \nabla^2 \rho(\mathbf{r}) \ln \frac{\rho(\mathbf{r})}{\rho_0(\mathbf{r})}\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
prorho (np.ndarray((N,), dtype=float), optional) – Electron density of promolecule on grid of N points, by default None.
rho_lapl (np.ndarray((N,), dtype=float), optional) – Electron density laplacian on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- G2(rho=None, prorho=None, rho_lapl=None, prorho_lapl=None, omega=None)[source]
G2 density defined as:
\[g_2(\mathbf{r}) = \rho(\mathbf{r}) \left[ \frac{\nabla^2 \rho(\mathbf{r})}{\rho(\mathbf{r})} -\frac{\nabla^2 \rho_0(\mathbf{r})}{\rho_0(\mathbf{r})} \right]\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
prorho (np.ndarray((N,), dtype=float), optional) – Electron density of promolecule on grid of N points, by default None.
rho_lapl (np.ndarray((N,), dtype=float), optional) – Electron density laplacian on grid of N points, by default None.
prorho_lapl (np.ndarray((N,), dtype=float), optional) – Electron density laplacian of promolecule on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- G3(rho=None, prorho=None, rho_grad=None, prorho_grad=None, omega=None)[source]
G3 density defined as:
\[g_3(\mathbf{r}) \equiv {}^r_F i(\mathbf{r}) = \rho(\mathbf{r}) \left\vert \frac{\nabla \rho(\mathbf{r})}{\rho(\mathbf{r})} -\frac{\nabla \rho_0(\mathbf{r})}{\rho_0(\mathbf{r})} \right\vert^2\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
prorho (np.ndarray((N,), dtype=float), optional) – Electron density of promolecule on grid of N points, by default None.
rho_grad (np.ndarray((N,), dtype=float), optional) – Electron density graident on grid of N points, by default None.
prorho_grad (np.ndarray((N,), dtype=float), optional) – Electron density graident of promolecule on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- GBP_entropy(rho=None, ts=None, tTF=None, k=1.0, omega=None)[source]
Ghosh-Berkowitz-Parr (GBP) entropy density \(s_{GBP}\) defined as:
\[s_{GBP} = \frac{3}{2}k\rho(\mathbf{r}) \left[ c+\ln \frac{t(\mathbf{r};\rho)}{t_{TF}(\mathbf{r};\rho)} \right]\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
ts (np.ndarray((N,), dtype=float), optional) – Single particle kenetic Electron density on grid of N points, by default None.
tTF (np.ndarray((N,), dtype=float), optional) – Thomas-Fermi kinetic energy density on grid of N points, by default None.
k (float, optional) – Boltzmann constant, by default 1.0 for convenience.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- alternative_fisher_information(rho=None, rho_lapl=None, omega=None)[source]
Alternative Fisher information density defined as:
\[I^{\prime}_F = \nabla^2 \rho(\mathbf{r}) \ln \rho(\mathbf{r})\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
rho_lapl (np.ndarray((N,), dtype=float), optional) – Electron density laplacian on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- fisher_information(rho=None, rho_grad_norm=None, omega=None)[source]
Fisher information density defined as:
\[i_F = \frac{|\nabla \rho(\mathbf{r})|^2}{\rho(\mathbf{r})}\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
rho_grad_norm (np.ndarray((N,), dtype=float), optional) – Electron density graident norm on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- relative_alternative_fisher_information(rho=None, prorho=None, rho_lapl=None, omega=None)[source]
Relative alternative Fisher information density defined as:
\[{}^r_F i^{\prime}(\mathbf{r}) = \nabla^2 \rho(\mathbf{r}) \ln \frac{\rho(\mathbf{r})}{\rho_0(\mathbf{r})}\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
prorho (np.ndarray((N,), dtype=float), optional) – Electron density of promolecule on grid of N points, by default None.
rho_lapl (np.ndarray((N,), dtype=float), optional) – Electron density laplacian on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- relative_fisher_information(rho=None, prorho=None, rho_grad=None, prorho_grad=None, omega=None)[source]
Relative Fisher information density defined as:
\[{}^r_F i(\mathbf{r}) = \rho(\mathbf{r}) \left\vert \frac{\nabla \rho(\mathbf{r})}{\rho(\mathbf{r})} -\frac{\nabla \rho_0(\mathbf{r})}{\rho_0(\mathbf{r})} \right\vert^2\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
prorho (np.ndarray((N,), dtype=float), optional) – Electron density of promolecule on grid of N points, by default None.
rho_grad (np.ndarray((N,), dtype=float), optional) – Electron density graident on grid of N points, by default None.
prorho_grad (np.ndarray((N,), dtype=float), optional) – Electron density graident of promolecule on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- relative_rho_power(n=2, rho=None, prorho=None, omega=None)[source]
Relative Renyi entropy density defined as:
\[r^r_n = \frac{\rho^n(\mathbf{r})}{\rho^{n-1}_0(\mathbf{r})}\]- Parameters:
n (int, optional) – Order of relative Renyi entropy, by default 2.
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
prorho (np.ndarray((N,), dtype=float), optional) – Electron density of promolecule on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- relative_shannon_entropy(rho=None, prorho=None, omega=None)[source]
Relative Shannon entropy density defined as:
\[s^r_S = -\rho(\mathbf{r})\ln \frac{\rho(\mathbf{r})}{\rho_0(\mathbf{r})}\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
prorho (np.ndarray((N,), dtype=float), optional) – Electron density of promolecule on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- rho_power(n=2, rho=None, omega=None)[source]
Electron density of power n defined as \(\rho(\mathbf{r})^n\).
- Parameters:
n (int, optional) – Order of rho power, by default 2.
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
- shannon_entropy(rho=None, omega=None)[source]
Shannon entropy density defined as:
\[s_S = -\rho(\mathbf{r}) \ln \rho(\mathbf{r})\]- Parameters:
rho (np.ndarray((N,), dtype=float), optional) – Electron density on grid of N points, by default None.
omega (np.ndarray((N,), dtype=float), optional) – Sharing function of single atom, by default None.
- Returns:
ita_density – Information theory density on grid of N points.
- Return type:
np.ndarray((N,), dtype=float)
3.6. pyscf.ita.ita module
Information-theoretic approach (ITA) module.
References
The following references are for the atoms-in-molecules module.
- class pyscf.ita.ita.ITA(method, grids, rung=[3, 3], spin='ab', promolecule=None)[source]
Bases:
object
Information-Theoretic Approch (ITA) class.
Examples:
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1.1') >>> mf = scf.HF(mol) >>> mf.kernel() >>> grids = dft.Grids(mol) >>> grids.build() >>> ita = ITA() >>> ita.method = mf >>> ita.grids = grids >>> ita.build()
- G1(grids_weights=None, ita_density=None)[source]
G1 density defined as:
\[G_1(\mathbf{r}) = \sum_A \int \nabla^2 \rho_A(\mathbf{r}) \ln \frac{\rho_A(\mathbf{r})}{\rho^0_A(\mathbf{r})} d\mathbf{r}\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- G2(grids_weights=None, ita_density=None)[source]
G2 density defined as:
\[G_2(\mathbf{r}) = \sum_A \int \rho_A(\mathbf{r}) \left[ \frac{\nabla^2 \rho_A(\mathbf{r})}{\rho_A(\mathbf{r})} -\frac{\nabla^2 \rho^0_A(\mathbf{r})}{\rho^0_A(\mathbf{r})} \right] d\mathbf{r}\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- G3(grids_weights=None, ita_density=None)[source]
G3 density defined as:
\[G_3(\mathbf{r}) = \sum_A \int \rho_A(\mathbf{r}) \left\vert \frac{\nabla \rho_A(\mathbf{r})}{\rho_A(\mathbf{r})} -\frac{\nabla \rho^0_A(\mathbf{r})}{\rho^0_A(\mathbf{r})} \right\vert^2 d\mathbf{r}\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- GBP_entropy(grids_weights=None, ita_density=None)[source]
Ghosh-Berkowitz-Parr (GBP) entropy \(S_{GBP}\) defined as:
\[S_{GBP} = \int \frac{3}{2}k\rho(\mathbf{r}) \left[ c+\ln \frac{t(\mathbf{r};\rho)}{t_{TF}(\mathbf{r};\rho)} \right] d\mathbf{r}\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- alternative_fisher_information(grids_weights=None, ita_density=None)[source]
Alternative Fisher information \(I^{\prime}_F\) defined as:
\[I^{\prime}_F = -\int \nabla^2 \rho(\mathbf{r}) \ln \rho(\mathbf{r}) d\mathbf{r}\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- batch_compute(ita_code=[], representation='electron density', partition='hirshfeld', filename='pyita.log')[source]
ITA batch calcuation.
- Parameters:
ita (ITA) – Instance of ITA class.
ita_code (List[int]) – List of ITA code to calculate.
representation (('electron density' | 'shape function' | 'atoms in molecules')) – Type of representation, by default ‘electron density’.
partition (('hirshfeld' | 'bader' | 'becke'), optional) – Atoms in molecule partition method, by default ‘hirshfeld’.
filename (str, optional) – File path and name of output, by default ‘pyita.log’
- build(method=None, grids=None, rung=None, spin=None, promolecule=None)[source]
Method to build the class.
- Parameters:
method (PyscfMethod, optional) – Pyscf scf method or post-scf method instance, by deafult None.
grids (Grids, optional) – Pyscf Grids instance, by deafult None.
rung (List[int, int], optional) – List of molecule and promolecule derivate+1 level, by default [3,0].
spin (('ab' | 'a' | 'b' | 'm'), optional) – Type of density to compute; either total, alpha-spin, beta-spin, or magnetization density, by default None.
promolecule (ProMolecule, optional) – ProMolecule instance, by default None.
- fisher_information(grids_weights=None, ita_density=None)[source]
Fisher information \(I_F\) defined as:
\[I_F = \int \frac{|\nabla \rho(\mathbf{r})|^2}{\rho(\mathbf{r})} d\mathbf{r}\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- onicescu_information(n=2, grids_weights=None, ita_density=None)[source]
Onicescu information \(E_n\) defined as:
\[E_n = \frac{1}{n-1} \int \rho(\mathbf{r})^n d\mathbf{r}\]- Parameters:
n (int, optional) – Order of the Onicescu information, by default 2.
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- relative_alternative_fisher_information(grids_weights=None, ita_density=None)[source]
Relative alternative Fisher information defined as:
\[{}^r_F I^{\prime}(\mathbf{r}) = \int \nabla^2 \rho(\mathbf{r}) \ln \frac{\rho(\mathbf{r})}{\rho_0(\mathbf{r})} dr\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- relative_fisher_information(grids_weights=None, ita_density=None)[source]
Relative Fisher information defined as:
\[{}^r_F I(\mathbf{r}) = \int \rho(\mathbf{r}) \left\vert \frac{\nabla \rho(\mathbf{r})}{\rho(\mathbf{r})} -\frac{\nabla \rho_0(\mathbf{r})}{\rho_0(\mathbf{r})} \right\vert^2 dr\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- relative_onicescu_information(n=2, grids_weights=None, ita_density=None)[source]
Relative Tsallis entropy \(E^r_n\) defined as:
\[E^r_n = \frac{1}{n-1} \int \frac{\rho^n(\mathbf{r})}{\rho^{n-1}_0(\mathbf{r})} d\mathbf{r}\]- Parameters:
n (int, optional) – Order of the relative Renyi entropy, by default 2.
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- relative_renyi_entropy(n=2, grids_weights=None, ita_density=None)[source]
Relative Renyi entropy \(R^r_n\) defined as:
\[R^r_n = \frac{1}{1-n} \ln \left[ \int \frac{\rho^n(\mathbf{r})}{\rho^{n-1}_0(\mathbf{r})} d\mathbf{r} \right]\]- Parameters:
n (int, optional) – Order of the relative Renyi entropy, by default 2.
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- relative_shannon_entropy(grids_weights=None, ita_density=None)[source]
Relative Shannon entropy \(S^r_S\) defined as:
\[S^r_S = \int \rho(\mathbf{r})\ln \frac{\rho(\mathbf{r})}{\rho_0(\mathbf{r})} d\mathbf{r}\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- relative_tsallis_entropy(n=2, grids_weights=None, ita_density=None)[source]
Relative Tsallis entropy \(T^r_n\) defined as:
\[T^r_n = \frac{1}{n-1} \left[ 1- \int \frac{\rho^n(\mathbf{r})}{\rho^{n-1}_0(\mathbf{r})} d\mathbf{r} \right]\]- Parameters:
n (int, optional) – Order of the relative Renyi entropy, by default 2.
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- renyi_entropy(n=2, grids_weights=None, ita_density=None)[source]
Renyi entropy \(R_n\) defined as:
\[R_n = \frac{1}{1-n} \ln \left[ \int \rho(\mathbf{r})^n d\mathbf{r} \right]\]- Parameters:
n (int, optional) – Order of the Renyi entropy, by default 2.
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- rho_power(n=2, grids_weights=None, ita_density=None)[source]
Electron density of power n defined as \(\int \rho(\mathbf{r})^n d\mathbf{r}\).
- Parameters:
n (int, optional) – Order of rho power, by default 2.
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- shannon_entropy(grids_weights=None, ita_density=None)[source]
Shannon entropy \(S_S\) defined as:
\[S_S = -\int \rho(\mathbf{r}) \ln \rho(\mathbf{r}) dr\]- Parameters:
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
- tsallis_entropy(n=2, grids_weights=None, ita_density=None)[source]
Tsallis entropy \(T_n\) defined as:
\[T_n = \frac{1}{n-1} \left[ 1- \int \rho(\mathbf{r})^n d\mathbf{r} \right]\]- Parameters:
n (int, optional) – Order of the Tsallis entropy, by default 2.
grids_weights (np.ndarray((N,), dtype=float), optional) – Grids weights on N points, by default None.
ita_density (np.ndarray((N,), dtype=float), optional) – ITA density, by default None.
- Returns:
result – Scalar ITA result.
- Return type:
float
3.7. pyscf.ita.script module
- pyscf.ita.script.batch_compute(ita, ita_code=[], representation='electron density', partition='hirshfeld', filename='ita.log')[source]
ITA batch calcuation.
- Parameters:
ita (ITA) – Instance of ITA class.
ita_code (List[int]) – List of ITA code to calculate.
representation (('electron density' | 'shape function' | 'atoms in molecules')) – Type of representation, by default ‘electron density’.
partition (('hirshfeld' | 'bader' | 'becke'), optional) – Atoms in molecule partition method, by default ‘hirshfeld’.
filename (str, optional) – File path and name of output, by default ‘ita.log’
- pyscf.ita.script.section_compute(ita, code, log, representation='electron density', partition='hirshfeld', **kwargs)[source]
Function to calculation a single ITA section.
- Parameters:
ita (ITA) – Instance of ITA class.
code (int) – Single ITA code.
log (Log) – Instance of Log class.
representation (('electron density' | 'shape function' | 'atoms in molecules')) – Type of representation, by default ‘electron density’.
partition (('hirshfeld' | 'bader' | 'becke'), optional) – Atoms in molecule partition method, by default ‘hirshfeld’.