Grider

class n2v.Grider[source]

Bases: object

Methods Summary

dft_grid_to_fock(value, Vpot)

For value on DFT spherical grid, Fock matrix is returned.

generate_dft_grid(Vpot)

Extracts DFT spherical grid and weights from wfn object

generate_grids(x, y, z)

Genrates Mesh from 3 separate linear spaces and flatten, needed for cubic grid.

get_basis_set_correction(grid)

grid_to_blocks(grid[, basis])

Generate list of blocks to allocate given grid

on_grid_ao(coeff[, grid, basis, Vpot])

Generates a quantity on the grid given its ao representation. *This is the most general function for basis to grid transformation.

on_grid_density([grid, Da, Db, Vpot])

Generates Density given grid

on_grid_esp([Da, Db, grid, Vpot, wfn])

Generates EXTERNAL/ESP/HARTREE and Fermi Amaldi Potential on given grid

on_grid_grad_phi([Ca, Cb, grid, Vpot])

Generates laplacian of molecular orbitals

on_grid_lap_phi([Ca, Cb, grid, Vpot])

Generates laplacian of molecular orbitals

on_grid_orbitals([Ca, Cb, grid, Vpot])

Generates orbitals given grid

on_grid_vxc([func_id, grid, Da, Db, Vpot])

Generates Vxc given grid

Methods Documentation

dft_grid_to_fock(value, Vpot)[source]

For value on DFT spherical grid, Fock matrix is returned. VFock_ij = int dx phi_i(x) phi_j(x) value(x)

value: np.ndarray of shape (npoint, ).

Vpot:psi4.core.VBase

Vpotential object with info about grid. Provides DFT spherical grid. Only comes to play if no grid is given.

VFock: np.ndarray of shape (nbasis, nbasis)

generate_dft_grid(Vpot)[source]

Extracts DFT spherical grid and weights from wfn object

Parameters

Vpot (psi4.core.VBase) – Vpot object with dft grid data

Returns

dft_grid – Numpy arrays corresponding to x,y,z, and w. Shape: (4, npoints)

Return type

list

generate_grids(x, y, z)[source]

Genrates Mesh from 3 separate linear spaces and flatten, needed for cubic grid.

Parameters

grid (tuple of three np.ndarray) – (x, y, z)

Returns

grid – shape (3, len(x)*len(y)*len(z)).

Return type

np.ndarray

get_basis_set_correction(grid)[source]
grid_to_blocks(grid, basis=None)[source]

Generate list of blocks to allocate given grid

Parameters
  • grid (np.ndarray) – Grid to be distributed into blocks Size: (3, npoints) for homogeneous grid

    (4, npoints) for inhomogenous grid to account for weights

  • basis (psi4.core.BasisSet; optional) – The basis set. If not given, it will use target wfn.basisset().

Returns

  • blocks (list) – List with psi4.core.BlockOPoints

  • npoints (int) – Total number of points (for one dimension)

  • points (psi4.core.{RKS, UKS}) – Points function to set matrices.

on_grid_ao(coeff, grid=None, basis=None, Vpot=None)[source]

Generates a quantity on the grid given its ao representation. *This is the most general function for basis to grid transformation.

Parameters
  • coeff (np.ndarray) – Vector/Matrix of quantity on ao basis. Shape: {(num_ao_basis, ), (num_ao_basis, num_ao_basis)}

  • grid (np.ndarray Shape: (3, npoints) or (4, npoints) or tuple for block_handler (return of grid_to_blocks)) – grid where density will be computed.

  • basis (psi4.core.BasisSet, optional) – The basis set. If not given it will use target wfn.basisset().

  • Vpot (psi4.core.VBase) – Vpotential object with info about grid. Provides DFT spherical grid. Only comes to play if no grid is given.

Returns

coeff_r – Quantity expressed by the coefficient on the given grid

Return type

np.ndarray Shape: (npoints, )

on_grid_density(grid=None, Da=None, Db=None, Vpot=None)[source]

Generates Density given grid

Parameters
  • Da, Db (np.ndarray) – Alpha, Beta densities. Shape: (num_ao_basis, num_ao_basis)

  • grid (np.ndarray Shape: (3, npoints) or (4, npoints) or tuple for block_handler (return of grid_to_blocks)) – grid where density will be computed.

  • Vpot (psi4.core.VBase) – Vpotential object with info about grid. Provides DFT spherical grid. Only comes to play if no grid is given.

Returns

density – Density on the given grid.

Return type

np.ndarray Shape: (ref, npoints)

on_grid_esp(Da=None, Db=None, grid=None, Vpot=None, wfn=None)[source]

Generates EXTERNAL/ESP/HARTREE and Fermi Amaldi Potential on given grid

Parameters
  • Da,Db (np.ndarray, opt, shape (nbf, nbf)) – The electron density in the denominator of Hartee potential. If None, the original density matrix will be used.

  • grid (np.ndarray Shape: (3, npoints) or (4, npoints) or tuple for block_handler (return of grid_to_blocks)) – grid where density will be computed.

  • Vpot (psi4.core.VBase) – Vpotential object with info about grid. Provides DFT spherical grid. Only comes to play if no grid is given.

Returns

vext, hartree, esp, v_fa – External, Hartree, ESP, and Fermi Amaldi potential on the given grid Shape: (npoints, )

Return type

np.ndarray

on_grid_grad_phi(Ca=None, Cb=None, grid=None, Vpot=None)[source]

Generates laplacian of molecular orbitals

Ca, Cb: np.ndarray

Alpha, Beta Orbital Coefficient Matrix. Shape: (num_ao_basis, num_ao_basis)

grid: np.ndarray Shape: (3, npoints) or (4, npoints) or tuple for block_handler (return of grid_to_blocks)

grid where density will be computed.

Vpot: psi4.core.VBase

Vpotential object with info about grid. Provides DFT spherical grid. Only comes to play if no grid is given.

grad_phi: List[np.ndarray]. Where array is of shape (npoints, ref)

Gradient of molecular orbitals on the grid

on_grid_lap_phi(Ca=None, Cb=None, grid=None, Vpot=None)[source]

Generates laplacian of molecular orbitals

Ca, Cb: np.ndarray

Alpha, Beta Orbital Coefficient Matrix. Shape: (num_ao_basis, num_ao_basis)

grid: np.ndarray Shape: (3, npoints) or (4, npoints) or tuple for block_handler (return of grid_to_blocks)

grid where density will be computed.

Vpot: psi4.core.VBase

Vpotential object with info about grid. Provides DFT spherical grid. Only comes to play if no grid is given.

lap_phi: List[np.ndarray]. Where array is of shape (npoints, ref)

Laplacian of molecular orbitals on the grid

on_grid_orbitals(Ca=None, Cb=None, grid=None, Vpot=None)[source]

Generates orbitals given grid

Parameters
  • Ca, Cb (np.ndarray) – Alpha, Beta Orbital Coefficient Matrix. Shape: (num_ao_basis, num_ao_basis)

  • grid (np.ndarray Shape: (3, npoints) or (4, npoints) or tuple for block_handler (return of grid_to_blocks)) – grid where density will be computed

  • Vpot (psi4.core.VBase) – Vpotential object with info about grid. Provides DFT spherical grid. Only comes to play if no grid is given.

Returns

orbitals – Orbitals on the given grid of size . Shape: (nbasis, npoints, ref)

Return type

np.ndarray

on_grid_vxc(func_id=1, grid=None, Da=None, Db=None, Vpot=None)[source]

Generates Vxc given grid

Parameters
  • Da, Db (np.ndarray) – Alpha, Beta densities. Shape: (num_ao_basis, num_ao_basis)

  • func_id (int) – Functional ID associated with Density Functional Approximationl. Full list of functionals: https://www.tddft.org/programs/libxc/functionals/

  • grid (np.ndarray Shape: (3, npoints) or (4, npoints) or tuple for block_handler (return of grid_to_blocks)) – grid where density will be computed.

  • Vpot (psi4.core.VBase) – Vpotential object with info about grid. Provides DFT spherical grid. Only comes to play if no grid is given.

Returns

VXC – Exchange correlation potential on the given grid Shape: (npoints, )

Return type

np.ndarray