gbasis.evals package#

Collection of modules that compute evaluations of different properties of the contractions.

Submodules#

gbasis.evals.density module#

Density Evaluation.

gbasis.evals.density.evaluate_density(one_density_matrix, basis, points, transform=None, threshold=1e-08)#

Return the density of the given basis set at the given points.

\[\rho(\mathbf{r}) = \sum_{ij} \gamma_{ij} \phi_i(\mathbf{r}) \phi_j(\mathbf{r})\]

where \(\mathbf{r}\) is the point at which the density is evaluated, \(\gamma_{ij}\) is the one-electron density matrix, and \(\phi_i\) is the \(i\)-th basis function.

Parameters#

one_density_matrixnp.ndarray(K_orbs, K_orbs)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

thresholdfloat, optional

The absolute value below which negative density values are acceptable. Any negative density value with an absolute value smaller than this threshold will be set to zero.

Returns#

densitynp.ndarray(N,)

Density evaluated at N grid points.

gbasis.evals.density.evaluate_density_gradient(one_density_matrix, basis, points, transform=None, deriv_type='general')#

Return the gradient of the density evaluated at the given points.

\[\begin{split}\nabla \rho(\mathbf{r}) = \begin{bmatrix} \frac{\partial}{\partial x} \rho(\mathbf{r})\\\\ \frac{\partial}{\partial y} \rho(\mathbf{r})\\\\ \frac{\partial}{\partial z} \rho(\mathbf{r}) \end{bmatrix} \end{split}\]

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

deriv_type“general” or “direct”

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

Returns#

density_gradientnp.ndarray(N, 3)

Gradient of the density evaluated at N grid points.

gbasis.evals.density.evaluate_density_hessian(one_density_matrix, basis, points, transform=None, deriv_type='general')#

Return the Hessian of the density evaluated at the given points.

\[\begin{split}H[\rho(\mathbf{r})] = \begin{bmatrix} \frac{\partial^2}{\partial x^2} \rho(\mathbf{r}) & \frac{\partial^2}{\partial x \partial y} \rho(\mathbf{r}) & \frac{\partial^2}{\partial x \partial z} \rho(\mathbf{r})\\\\ \frac{\partial^2}{\partial y \partial x} \rho(\mathbf{r}) & \frac{\partial^2}{\partial y^2} \rho(\mathbf{r})& \frac{\partial^2}{\partial y \partial z} \rho(\mathbf{r})\\\\ \frac{\partial^2}{\partial z \partial x} \rho(\mathbf{r}) & \frac{\partial^2}{\partial z \partial y} \rho(\mathbf{r})& \frac{\partial^2}{\partial z^2} \rho(\mathbf{r})\\ \end{bmatrix}\end{split}\]

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

deriv_type“general” or “direct”

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

Returns#

density_hessiannp.ndarray(N, 3, 3)

Hessian of the density evaluated at N grid points. Dimension 0 corresponds to the point, ordered as in points. Dimensions 1, 2 correspond to the dimensions (x, y, z) in which the derivative of density was calculated.

gbasis.evals.density.evaluate_density_laplacian(one_density_matrix, basis, points, transform=None, deriv_type='general')#

Return the Laplacian of the density evaluated at the given points.

\[\nabla^2 \rho(\mathbf{r}) = \frac{\partial^2}{\partial x^2} \rho(\mathbf{r}) + \frac{\partial^2}{\partial y^2} \rho(\mathbf{r}) + \frac{\partial^2}{\partial z^2} \rho(\mathbf{r})\]

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

deriv_type“general” or “direct”

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

Returns#

density_laplaciannp.ndarray(N)

Laplacian of the density evaluated at N grid points.

gbasis.evals.density.evaluate_density_using_evaluated_orbs(one_density_matrix, orb_eval)#

Return the evaluation of the density given the evaluated orbitals.

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given orbitals.

orb_evalnp.ndarray(K_orb, N)

Orbitals evaluated at \(N\) grid points. The set of orbitals must be the same basis set used to build the one-electron density matrix.

Returns#

densitynp.ndarray(N,)

Density evaluated at N grid points.

Raises#

TypeError

If orb_eval is not a 2-dimensional numpy array with dtype float. If one_density_matrix is not a 2-dimensional numpy array with dtype float.

ValueError

If one_density_matrix is not square. If the number of columns (or rows) of one_density_matrix is not equal to the number of rows of the orbital evaluations.

gbasis.evals.density.evaluate_deriv_density(orders, one_density_matrix, basis, points, transform=None, deriv_type='general')#

Return the derivative of density of the given transformed basis set at the given points.

\[\begin{split}&\frac{\partial^{L_x + L_y + L_z}}{\partial x^{L_x} \partial y^{L_y} \partial z^{L_z}} \rho(\mathbf{r})\\\\ &= \sum_{l_x=0}^{L_x} \sum_{l_y=0}^{L_y} \sum_{l_z=0}^{L_z} \binom{L_x}{l_x} \binom{L_y}{l_y} \binom{L_z}{l_z} \sum_{ij} \gamma_{ij} \frac{\partial^{l_x + l_y + l_z} \rho(\mathbf{r})}{\partial x^{l_x} \partial y^{l_y} \partial z^{l_z}} \frac{ \partial^{L_x + L_y + L_z - l_x - l_y - l_z} \rho(\mathbf{r}) }{ \partial x^{L_x - l_x} \partial y^{L_y - l_y} \partial z^{L_z - l_z} }\end{split}\]

where \(L_x, L_y, L_z\) are the orders of the derivative relative to the \(x, y, \text{and} z\) components, respectively.

Parameters#

ordersnp.ndarray(3,)

Orders of the derivative.

one_density_matrixnp.ndarray(K_orbs, K_orbs)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

deriv_type“general” or “direct”

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

Returns#

density_derivnp.ndarray(N,)

Derivative of the density evaluated at N grid points.

gbasis.evals.density.evaluate_deriv_reduced_density_matrix(orders_one, orders_two, one_density_matrix, basis, points, transform=None, deriv_type='general')#

Return the derivative of the first-order reduced density matrix at the given points.

\[\begin{split}&\left. \frac{\partial^{p_x + p_y + p_z}}{\partial x_1^{p_x} \partial y_1^{p_y} \partial z_1^{p_z}} \frac{\partial^{q_x + q_y + q_z}}{\partial x_2^{q_x} \partial y_2^{q_y} \partial z_2^{q_z}} \gamma(\mathbf{r}_1, \mathbf{r}_2) \right|_{\mathbf{r}_1 = \mathbf{r}_2 = \mathbf{r}}\\\\ &= \sum_{ij} \gamma_{ij} \left. \frac{\partial^{p_x + p_y + p_z}}{\partial x_1^{p_x} \partial y_1^{p_y} \partial z_1^{p_z}} \phi_i(\mathbf{r}_1) \right|_{\mathbf{r}_1 = \mathbf{r}} \left. \frac{\partial^{q_x + q_y + q_z}}{\partial x_2^{q_x} \partial y_2^{q_y} \partial z_2^{q_z}} \phi_j(\mathbf{r}_2) \right|_{\mathbf{r}_2 = \mathbf{r}}\end{split}\]

where \(\mathbf{r}_1\) is the first point, \(\mathbf{r}_2\) is the second point, and \(\mathbf{r}\) is the point at which the derivative is evaluated.

Parameters#

orders_onenp.ndarray(3,)

Orders of the derivative for the first point, \(mathbf{r}_1\).

orders_twonp.ndarray(3,)

Orders of the derivative for the second point, \(mathbf{r}_1\).

one_density_matrixnp.ndarray(K_orbs, K_orbs)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

deriv_type“general” or “direct”

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

Returns#

deriv_reduced_density_matrixnp.ndarray(N,)

Derivative of the first-order reduced density matrix evaluated at N grid points.

gbasis.evals.density.evaluate_general_kinetic_energy_density(one_density_matrix, basis, points, alpha, transform=None, deriv_type='general')#

Return evaluations of general form of the kinetic energy density at the given points.

\[t_{\alpha} (\mathbf{r}) = t_+(\mathbf{r}) + \alpha \nabla^2 \rho(\mathbf{r})\]

where \(t_+\) is the positive definite kinetic energy density.

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

alphafloat

Parameter of the general form of the kinetic energy density.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

deriv_type“general” or “direct”

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

Returns#

general_kinetic_energy_densitynp.ndarray(N,)

General kinetic energy density of the given transformed basis set evaluated at N grid points.

Raises#

TypeError

If alpha is not an integer or a float.

gbasis.evals.density.evaluate_posdef_kinetic_energy_density(one_density_matrix, basis, points, transform=None, deriv_type='general', threshold=1e-08)#

Return evaluations of positive definite kinetic energy density at the given points.

\[\begin{split}\begin{split} t_+ (\mathbf{r}) &= \frac{1}{2} \left. \nabla_{\mathbf{r}_1} \cdot \nabla_{\mathbf{r}_2} \gamma(\mathbf{r}_1, \mathbf{r}_2) \right|_{\mathbf{r}_1 = \mathbf{r}_2 = \mathbf{r}}\\ &= \frac{1}{2} \left( \frac{\partial^2}{\partial x_1 \partial x_2} \gamma(\mathbf{r}_1, \mathbf{r}_2) + \frac{\partial^2}{\partial y_1 \partial y_2} \gamma(\mathbf{r}_1, \mathbf{r}_2) + \frac{\partial^2}{\partial z_1 \partial z_2} \gamma(\mathbf{r}_1, \mathbf{r}_2) \right)_{\mathbf{r}_1 = \mathbf{r}_2 = \mathbf{r}}\\ \end{split}\end{split}\]

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont), optional

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions.

deriv_type“general” or “direct”, optional

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

thresholdfloat, optional

The absolute value below which negative density values are acceptable. Any negative density value with an absolute value smaller than this threshold will be set to zero.

Returns#

posdef_kindetic_energy_densitynp.ndarray(N,)

Positive definite kinetic energy density of the given transformed basis set evaluated at N grid points.

gbasis.evals.electrostatic_potential module#

Module for computing electrostatic potential integrals.

gbasis.evals.electrostatic_potential.electrostatic_potential(basis, one_density_matrix, points, nuclear_coords, nuclear_charges, transform=None, threshold_dist=0.0)#

Return the electrostatic potentials of the basis set in the Cartesian form.

\[- \left( - \sum_A \frac{Z_A}{|\mathbf{R}_C - \mathbf{R}_A|} + \sum_{ab} \gamma_{ab} \int \phi_a(\mathbf{r}) \frac{-1}{|\mathbf{r} - \mathbf{R}_C|} \phi_b(\mathbf{r}) d\mathbf{r} \right)\]

where \(\mathbf{R}_C\) is the coordinate of a unitary point charge, \(\mathbf{R}_A\) is the coordinate of the nucleus \(A\), \(Z_A\) its charge, and \(\gamma_{ab}\) is the one-electron density matrix.

Parameters#

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

one_density_matrixnp.ndarray(K_orbs, K_orbs)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

nuclear_coordsnp.ndarray(N_nuc, 3)

Cartesian coordinates of each atom. Rows correspond to the atoms and columns correspond to the \(x, y, \text{and} z\) components.

nuclear_chargesnp.ndarray(N_nuc)

Charges of each atom.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

threshold_dist{float, 0.0}

Threshold for rejecting nuclei whose distances to the points are less than the provided value. i.e. nuclei that are closer to the point than the threshold are discarded when computing the electrostatic potential of the point. Default value is 0.0, i.e. no nuclei are discarded.

Returns#

arraynp.ndarray(N)

Electrostatic potential evaluated at the given points.

Raises#

TypeError

If one_density_matrix is not a two-dimensional numpy array. If nuclear_coords is not a two-dimensional numpy array with 3 columns. If nuclear_charges is not a one-dimensional numpy array. If threshold_dist is not an integer or float.

ValueError

If one_density_matrix is not a symmetric (square) matrix. If number of rows in nuclear_coords is not equal to the number of elements in nuclear_charges. If threshold_dist is less than 0.

gbasis.evals.eval module#

Functions for evaluating Gaussian contractions.

class gbasis.evals.eval.Eval(contractions)#

Bases: BaseOneIndex

Class for evaluating Gaussian contractions and their linear combinations.

Dimension 0 of the returned array is associated with a contracted Gaussian (or a linear combination of a set of Gaussians).

Attributes#

_axes_contractionstuple of tuple of GeneralizedContractionShell

Contractions that are associated with each index of the array. Each tuple of GeneralizedContractionShell corresponds to an index of the array.

contractionstuple of GeneralizedContractionShell

Contractions that are associated with the first index of the array. Property of Eval.

Methods#

__init__(self, contractions)

Initialize.

construct_array_contraction(contraction, points)np.ndarray(M, L_cart, N)

Return the evaluations of the given Cartesian contractions at the given coordinates. M is the number of segmented contractions with the same exponents (and angular momentum). L_cart is the number of Cartesian contractions for the given angular momentum. N is the number of coordinates at which the contractions are evaluated.

construct_array_cartesian(self, points)np.ndarray(K_cart, N)

Return the evaluations of the Cartesian contractions of the instance at the given coordinates. K_cart is the total number of Cartesian contractions within the instance. N is the number of coordinates at which the contractions are evaluated.

construct_array_spherical(self, points)np.ndarray(K_sph, N)

Return the evaluations of the spherical contractions of the instance at the given coordinates. K_sph is the total number of spherical contractions within the instance. N is the number of coordinates at which the contractions are evaluated.

construct_array_mix(self, coord_types, points)np.ndarray(K_cont, N)

Return the evaluatations of the contraction in the given coordinate system. K_cont is the total number of contractions within the given basis set. N is the number of coordinates at which the contractions are evaluated.

construct_array_lincomb(self, transform, coord_type, points)np.ndarray(K_orbs, N)

Return the evaluations of the linear combinations of contractions in the given coordinate system. K_orbs is the number of basis functions produced after the linear combinations. N is the number of coordinates at which the contractions are evaluated.

static construct_array_contraction(contractions, points)#

Return the evaluations of the given contractions at the given coordinates.

Parameters#

contractionsGeneralizedContractionShell

Contracted Cartesian Gaussians (of the same shell) that will be used to construct an array.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

Returns#

array_contractionnp.ndarray(M, L_cart, N)

Evaluations of the given Cartesian contractions at the given points. Dimension 0 corresponds to segmented contractions within the given generalized contraction (same exponents and angular momentum, but different coefficients). M is the number of segmented contractions with the same exponents (and angular momentum). Dimension 1 corresponds to angular momentum vector. L_cart is the number of Cartesian contractions for the given angular momentum. Dimension 2 corresponds to coordinates at which the contractions are evaluated. N is the number of coordinates at which the contractions are evaluated.

Raises#

TypeError

If contractions is not a GeneralizedContractionShell instance. If points is not a two-dimensional numpy array with 3 columns.

Note#

Since all of the keyword arguments of construct_array_cartesian, construct_array_spherical, and construct_array_lincomb are ultimately passed down to this method, all of the mentioned methods must be called with the keyword arguments points and orders.

gbasis.evals.eval.evaluate_basis(basis, points, transform=None)#

Evaluate the basis set in the given coordinate system at the given points.

Parameters#

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

Returns#

eval_arraynp.ndarray(K, N)

Evaluations of the basis functions at the given points. If keyword argument transform is provided, then the transformed basis functions will be evaluated at the given points. K is the total number of basis functions within the given basis set. N is the number of coordinates at which the contractions are evaluated.

gbasis.evals.eval_deriv module#

Functions for evaluating Gaussian primitives.

class gbasis.evals.eval_deriv.EvalDeriv(contractions)#

Bases: BaseOneIndex

Class for evaluating Gaussian contractions and their linear combinations.

Dimension 0 of the returned array is associated with a contracted Gaussian (or a linear combination of a set of Gaussians).

Attributes#

_axes_contractionstuple of tuple of GeneralizedContractionShell

Contractions that are associated with each index of the array. Each tuple of GeneralizedContractionShell corresponds to an index of the array.

contractionstuple of GeneralizedContractionShell

Contractions that are associated with the first index of the array. Property of EvalDeriv.

Methods#

__init__(self, contractions)

Initialize.

construct_array_contraction(contraction, points, orders)np.ndarray(M, L_cart, N)

Return the evaluations of the given Cartesian contractions at the given coordinates. M is the number of segmented contractions with the same exponents (and angular momentum). L_cart is the number of Cartesian contractions for the given angular momentum. N is the number of coordinates at which the contractions are evaluated.

construct_array_cartesian(self, points, orders)np.ndarray(K_cart, N)

Return the evaluations of the derivatives of the Cartesian contractions of the instance at the given coordinates. K_cart is the total number of Cartesian contractions within the instance. N is the number of coordinates at which the contractions are evaluated.

construct_array_spherical(self, points, orders)np.ndarray(K_sph, N)

Return the evaluations of the derivatives of the spherical contractions of the instance at the given coordinates. K_sph is the total number of spherical contractions within the instance. N is the number of coordinates at which the contractions are evaluated.

construct_array_mix(self, coord_types, points, orders)np.ndarray(K_cont, N)

Return the array associated with all of the contraction in the given coordinate system. K_cont is the total number of contractions within the given basis set. N is the number of coordinates at which the contractions are evaluated.

construct_array_lincomb(self, transform, coord_type, points, orders)np.ndarray(K_orbs, N)

Return the evaluation of derivatives of the linear combinations of contractions in the given coordinate system. K_orbs is the number of basis functions produced after the linear combinations. N is the number of coordinates at which the contractions are evaluated.

static construct_array_contraction(contractions, points, orders, deriv_type='general')#

Return the array associated with a set of contracted Cartesian Gaussians.

Parameters#

contractionsGeneralizedContractionShell

Contracted Cartesian Gaussians (of the same shell) that will be used to construct an array.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

ordersnp.ndarray(3,)

Orders of the derivative.

deriv_type“general” or “direct”

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

Returns#

array_contractionnp.ndarray(M, L_cart, N)

Array associated with the given instance(s) of GeneralizedContractionShell. Dimension 0 corresponds to segmented contractions within the given generalized contraction (same exponents and angular momentum, but different coefficients). M is the number of segmented contractions with the same exponents (and angular momentum). Dimension 1 corresponds to angular momentum vector. L_cart is the number of Cartesian contractions for the given angular momentum. Dimension 2 corresponds to coordinates at which the contractions are evaluated. N is the number of coordinates at which the contractions are evaluated.

Raises#

TypeError

If contractions is not a GeneralizedContractionShell instance. If points is not a two-dimensional numpy array with 3 columns. If orders is not a one-dimensional numpy array with 3 elements.

ValueError

If orders has any negative numbers. If orders does not have dtype int.

Note#

Since all of the keyword arguments of construct_array_cartesian, construct_array_spherical, and construct_array_lincomb are ultimately passed down to this method, all of the mentioned methods must be called with the keyword arguments points and orders.

gbasis.evals.eval_deriv.evaluate_deriv_basis(basis, points, orders, transform=None, deriv_type='general')#

Evaluate the derivative of the basis set in the given coordinate system at the given points.

The derivative (to arbitrary orders) of a basis function is given by:

\[\frac{\partial^{m_x + m_y + m_z}}{\partial x^{m_x} \partial y^{m_y} \partial z^{m_z}} \phi (\mathbf{r})\]

where \(m_x, m_y, m_z\) are the orders of the derivative with respect to x, y, and z, \(\phi\) is the basis function (a generalized contraction shell), and \(\mathbf{r}_n\) are the coordinate of the points in space where the basis function is evaluated.

Parameters#

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

ordersnp.ndarray(3,)

Orders of the derivative. First element corresponds to the order of the derivative with respect to x. Second element corresponds to the order of the derivative with respect to y. Thirds element corresponds to the order of the derivative with respect to z.

transformnp.ndarray(K, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

deriv_type“general” or “direct”

Specification of derivative of contraction function in _deriv.py. “general” makes reference to general implementation of any order derivative function (_eval_deriv_contractions()) and “direct” makes reference to specific implementation of first and second order derivatives for generalized contraction (_eval_first_second_order_deriv_contractions()).

Returns#

eval_arraynp.ndarray(K, N)

Evaluations of the derivative of the basis functions at the given points. If keyword argument transform is provided, then the transformed basis functions will be evaluated at the given points. K is the total number of basis functions within the given basis set. N is the number of coordinates at which the contractions are evaluated.

gbasis.evals.stress_tensor module#

Module for computing properties related to the stress tensor.

gbasis.evals.stress_tensor.evaluate_ehrenfest_force(one_density_matrix, basis, points, alpha=1, beta=0, transform=None)#

Return the Ehrenfest force.

Ehrenfest force is the negative of the divergence of the stress tensor:

\[\begin{split}F_{j}(\mathbf{r} | \alpha, \beta) =&- \sum_i \frac{\partial}{\partial r_i} \boldsymbol{\sigma}_{ij}\\ =& \alpha \sum_i \left. \frac{\partial^3}{\partial r^2_i \partial r'_j} \gamma(\mathbf{r}, \mathbf{r}') \right|_{\mathbf{r} = \mathbf{r}'}\\ &- (1 - \alpha) \sum_i \left. \frac{\partial^3}{\partial r^2_i \partial r_j} \gamma(\mathbf{r}, \mathbf{r}) \right|_{\mathbf{r} = \mathbf{r}'} - (1 - 2\alpha) \sum_i \left. \frac{\partial^3}{\partial r_i \partial r_j \partial r'_i} \gamma(\mathbf{r}, \mathbf{r}) \right|_{\mathbf{r} = \mathbf{r}'}\\ &+ \frac{1}{2} \beta \left( \frac{\partial^3}{\partial r_j \partial x^2} + \frac{\partial^3}{\partial r_j \partial y^2} + \frac{\partial^3}{\partial r_j \partial z^2} \right) \rho(\mathbf{r})\end{split}\]

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

alpha{int, float}

First parameter of the stress tensor. Default value is 1.

beta{int, float}

Second parameter of the stress tensor. Default value is 0.

Returns#

ehrenfest_forcenp.ndarray(N, 3)

Ehrenfest force of the given density matrix evaluated at the given coordinates.

Raises#

TypeError

If alpha is not an integer or float. If beta is not an integer or float.

gbasis.evals.stress_tensor.evaluate_ehrenfest_hessian(one_density_matrix, basis, points, alpha=1, beta=0, transform=None, symmetric=False)#

Return the Ehrenfest Hessian.

Ehrenfest Hessian is the gradient of the Ehrenfest force:

\[\begin{split}H_{jk}(\mathbf{r} | \alpha, \beta) =& - \frac{\partial}{\partial r_k} F_j(\mathbf{r} | \alpha, \beta)\\ =& \alpha \sum_i \left( \frac{\partial^4}{\partial r^2_i \partial r_k \partial r'_j} \gamma(\mathbf{r}, \mathbf{r}') +\frac{\partial^4}{\partial r^2_i \partial r'_j \partial r'_k} \gamma(\mathbf{r}, \mathbf{r}') \right)_{\mathbf{r} = \mathbf{r}'}\\ &- (1 - \alpha) \sum_i \left( \frac{\partial^4}{\partial r^2_i \partial r_j \partial r_k} \gamma(\mathbf{r}, \mathbf{r}) + \frac{\partial^4}{\partial r^2_i \partial r_j \partial r'_k} \gamma(\mathbf{r}, \mathbf{r}) \right)_{\mathbf{r} = \mathbf{r}'}\\ &- (1 - 2\alpha) \sum_i \left( \frac{\partial^4}{\partial r_i \partial r_j \partial r_k \partial r'_i} \gamma(\mathbf{r}, \mathbf{r}) + \frac{\partial^4}{\partial r_i \partial r_j \partial r'_i \partial r'_k} \gamma(\mathbf{r}, \mathbf{r}) \right)_{\mathbf{r} = \mathbf{r}'}\\ &+ \frac{1}{2} \beta \left( \frac{\partial^4}{\partial r_j \partial r_k \partial x^2} + \frac{\partial^4}{\partial r_j \partial r_k \partial y^2} + \frac{\partial^4}{\partial r_j \partial r_k \partial z^2} \right) \rho(\mathbf{r})\end{split}\]

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

alpha{int, float}

First parameter of the stress tensor. Default value is 1.

beta{int, float}

Second parameter of the stress tensor. Default value is 0.

symmetric{True, False}

Flag for symmetrizing the Hessian. If True, then the Hessian is symmetrized by averaging it with its transpose. Default value is False.

Returns#

ehrenfest_hessiannp.ndarray(N, 3, 3)

Ehrenfest Hessian of the given density matrix evaluated at the given coordinates. Hessian is symmetrized if symmetric is True.

Raises#

TypeError

If alpha is not an integer or float. If beta is not an integer or float.

gbasis.evals.stress_tensor.evaluate_stress_tensor(one_density_matrix, basis, points, alpha=1, beta=0, transform=None)#

Return the stress tensor evaluated at the given coordinates.

Stress tensor is defined here as:

\[\begin{split}\boldsymbol{\sigma}_{ij}(\mathbf{r} | \alpha, \beta) =& -\frac{1}{2} \alpha \left( \frac{\partial^2}{\partial r_i \partial r'_j} \gamma(\mathbf{r}, \mathbf{r}') + \frac{\partial^2}{\partial r_j \partial r'_i} \gamma(\mathbf{r}, \mathbf{r}') \right)_{\mathbf{r} = \mathbf{r}'}\\ & +\frac{1}{2} (1 - \alpha) \left( \frac{\partial^2}{\partial r_i \partial r_j} \gamma(\mathbf{r}, \mathbf{r}) + \frac{\partial^2}{\partial r'_i \partial r'_j} \gamma(\mathbf{r}, \mathbf{r}') \right)_{\mathbf{r} = \mathbf{r}'}\\ & - \frac{1}{2} \delta_{ij} \beta \left. \nabla^2 \rho(\mathbf{r}) \right)\\ =& - \alpha \left. \frac{\partial^2}{\partial r_i \partial r'_j} \gamma(\mathbf{r}, \mathbf{r}') \right|_{\mathbf{r} = \mathbf{r}'} + (1 - \alpha) \left. \frac{\partial^2}{\partial r_i \partial r_j} \gamma(\mathbf{r}, \mathbf{r}) \right|_{\mathbf{r} = \mathbf{r}'} - \frac{1}{2} \delta_{ij} \beta \left. \nabla^2 \rho(\mathbf{r}) \right)\end{split}\]

Parameters#

one_density_matrixnp.ndarray(K_orb, K_orb)

One-electron density matrix in terms of the given basis set. If the basis is transformed using transform keyword, then the density matrix is assumed to be expressed with respect to the transformed basis set.

basislist/tuple of GeneralizedContractionShell

Shells of generalized contractions.

pointsnp.ndarray(N, 3)

Cartesian coordinates of the points in space (in atomic units) where the basis functions are evaluated. Rows correspond to the points and columns correspond to the \(x, y, \text{and} z\) components.

transformnp.ndarray(K_orbs, K_cont)

Transformation matrix from the basis set in the given coordinate system (e.g. AO) to linear combinations of contractions (e.g. MO). Transformation is applied to the left, i.e. the sum is over the index 1 of transform and index 0 of the array for contractions. Default is no transformation.

alpha{int, float}

First parameter of the stress tensor. Default value is 1.

beta{int, float}

Second parameter of the stress tensor. Default value is 0.

Returns#

stress_tensornp.ndarray(N, 3, 3)

Stress tensor of the given density matrix evaluated at the given points.

Raises#

TypeError

If alpha is not an integer or float. If beta is not an integer or float.