pyscf.pbc.tools package

Submodules

pyscf.pbc.tools.k2gamma module

Convert the k-sampled MO/integrals to corresponding Gamma-point supercell MO/integrals.

Zhihao Cui <zcui@caltech.edu> Qiming Sun <osirpt.sun@gmail.com>

See also the original implementation at https://github.com/zhcui/local-orbital-and-cdft/blob/master/k2gamma.py

pyscf.pbc.tools.k2gamma.double_translation_indices(kmesh)[source]

Indices to utilize the translation symmetry in the 2D matrix.

D[M,N] = D[N-M]

The return index maps the 2D subscripts to 1D subscripts.

D2 = D1[double_translation_indices()]

D1 holds all the symmetry unique elements in D2

pyscf.pbc.tools.k2gamma.get_phase(cell, kpts, kmesh=None, wrap_around=False)[source]

The unitary transformation that transforms the supercell basis k-mesh adapted basis.

pyscf.pbc.tools.k2gamma.k2gamma(kmf, kmesh=None)[source]

convert the k-sampled mean-field object to the corresponding supercell gamma-point mean-field object.

math:

C_{nu ‘ n’} = C_{vecRmu, veck m} = frac{1}{sqrt{N_{UC}}} e^{ii veckcdotvecR} C^{veck}_{mu m}

pyscf.pbc.tools.k2gamma.kpts_to_kmesh(cell, kpts, precision=None)[source]

Find the minimal k-points mesh to include all input kpts

pyscf.pbc.tools.k2gamma.mo_k2gamma(cell, mo_energy, mo_coeff, kpts, kmesh=None)[source]
pyscf.pbc.tools.k2gamma.to_supercell_ao_integrals(cell, kpts, ao_ints)[source]

Transform from the unitcell k-point AO integrals to the supercell gamma-point AO integrals.

pyscf.pbc.tools.k2gamma.to_supercell_mo_integrals(kmf, mo_ints)[source]

Transform from the unitcell k-point MO integrals to the supercell gamma-point MO integrals.

pyscf.pbc.tools.k2gamma.translation_map(nk)[source]

Generate [0 1 .. n ] [n 0 .. n-1] [n-1 n .. n-2] [… … …] [1 2 .. 0 ]

pyscf.pbc.tools.k2gamma.translation_vectors_for_kmesh(cell, kmesh, wrap_around=False)[source]

Translation vectors to construct super-cell of which the gamma point is identical to the k-point mesh of primitive cell

pyscf.pbc.tools.lattice module

pyscf.pbc.tools.lattice.get_ase_alkali_halide(A='Li', B='Cl')[source]
pyscf.pbc.tools.lattice.get_ase_atom(formula)[source]
pyscf.pbc.tools.lattice.get_ase_diamond_cubic(atom='C')[source]

Get the ASE atoms for cubic (8-atom) diamond unit cell.

pyscf.pbc.tools.lattice.get_ase_diamond_primitive(atom='C')[source]

Get the ASE atoms for primitive (2-atom) diamond unit cell.

pyscf.pbc.tools.lattice.get_ase_graphene(vacuum=5.0)[source]

Get the ASE atoms for primitive (2-atom) graphene unit cell.

pyscf.pbc.tools.lattice.get_ase_graphene_xxx(vacuum=5.0)[source]

Get the ASE atoms for primitive (2-atom) graphene unit cell.

pyscf.pbc.tools.lattice.get_ase_rocksalt(A='Li', B='Cl')[source]
pyscf.pbc.tools.lattice.get_ase_wurtzite(A='Zn', B='O')[source]
pyscf.pbc.tools.lattice.get_ase_zincblende(A='Ga', B='As')[source]
pyscf.pbc.tools.lattice.get_bandpath_fcc(ase_atom, npoints=30)[source]

pyscf.pbc.tools.make_test_cell module

pyscf.pbc.tools.make_test_cell.make_cell(L, mesh)[source]
pyscf.pbc.tools.make_test_cell.test_cell_cu_metallic(mesh=[9, 9, 9])[source]

Copper unit cell w/ special basis giving non-equal number of occupied orbitals per k-point

pyscf.pbc.tools.make_test_cell.test_cell_n0(L=5, mesh=[9, 9, 9])[source]
pyscf.pbc.tools.make_test_cell.test_cell_n1(L=5, mesh=[9, 9, 9])[source]
pyscf.pbc.tools.make_test_cell.test_cell_n2(L=5, mesh=[9, 9, 9])[source]
pyscf.pbc.tools.make_test_cell.test_cell_n3(mesh=[9, 9, 9])[source]

Take ASE Diamond structure, input into PySCF and run

pyscf.pbc.tools.make_test_cell.test_cell_n3_diffuse()[source]

Take ASE Diamond structure, input into PySCF and run

pyscf.pbc.tools.pbc module

pyscf.pbc.tools.pbc.cell_plus_imgs(cell, nimgs)[source]

Create a supercell via nimgs[i] in each +/- direction, as in get_lattice_Ls(). Note this function differs from super_cell() that super_cell only stacks the images in + direction.

Args:

cell : instance of Cell nimgs : (3,) array

Returns:

supcell : instance of Cell

pyscf.pbc.tools.pbc.cutoff_to_gs(a, cutoff)[source]

Deprecated. Replaced by function cutoff_to_mesh.

pyscf.pbc.tools.pbc.cutoff_to_mesh(a, cutoff)[source]

Convert KE cutoff to FFT-mesh

uses KE = k^2 / 2, where k_max ~ pi / grid_spacing

Args:
a(3,3) ndarray

The real-space cell lattice vectors. Each row represents a lattice vector.

cutofffloat

KE energy cutoff in a.u.

Returns:

mesh : (3,) array

pyscf.pbc.tools.pbc.fft(f, mesh)[source]

Perform the 3D FFT from real (R) to reciprocal (G) space.

After FFT, (u, v, w) -> (j, k, l). (jkl) is in the index order of Gv.

FFT normalization factor is 1., as in MH and in numpy.fft.

Args:
f(nx*ny*nz,) ndarray

The function to be FFT’d, flattened to a 1D array corresponding to the index order of cartesian_prod().

mesh(3,) ndarray of ints (= nx,ny,nz)

The number G-vectors along each direction.

Returns:
(nx*ny*nz,) ndarray

The FFT 1D array in same index order as Gv (natural order of numpy.fft).

pyscf.pbc.tools.pbc.fftk(f, mesh, expmikr)[source]

Perform the 3D FFT of a real-space function which is (periodic*e^{ikr}).

fk(k+G) = sum_r fk(r) e^{-i(k+G)r} = sum_r [f(k)e^{-ikr}] e^{-iGr}

pyscf.pbc.tools.pbc.get_coulG(cell, k=array([0., 0., 0.]), exx=False, mf=None, mesh=None, Gv=None, wrap_around=True, omega=None, **kwargs)[source]

Calculate the Coulomb kernel for all G-vectors, handling G=0 and exchange.

Args:
k(3,) ndarray

k-point

exxbool or str

Whether this is an exchange matrix element

mf : instance of SCF

Returns:
coulG(ngrids,) ndarray

The Coulomb kernel.

mesh(3,) ndarray of ints (= nx,ny,nz)

The number G-vectors along each direction.

omegafloat

Enable Coulomb kernel erf(|omega|*r12)/r12 if omega > 0 and erfc(|omega|*r12)/r12 if omega < 0. Note this parameter is slightly different to setting cell.omega for the treatment of exxdiv (at G0). cell.omega affects Ewald probe charge at G0. It is used mostly by RSH functionals for the long-range part of HF exchange. This parameter is used by range-separated JK builder and range-separated DF (and other range-separated integral methods) which require Ewald probe charge to be computed with regular Coulomb interaction (1/r12).

pyscf.pbc.tools.pbc.get_lattice_Ls(cell, nimgs=None, rcut=None, dimension=None, discard=True)[source]

Get the (Cartesian, unitful) lattice translation vectors for nearby images. The translation vectors can be used for the lattice summation.

Kwargs:
discard:

Drop less important Ls based on AO values on grid

pyscf.pbc.tools.pbc.get_monkhorst_pack_size(cell, kpts, tol=1e-05)[source]
pyscf.pbc.tools.pbc.gs_to_cutoff(a, gs)[source]

Deprecated. Replaced by function mesh_to_cutoff.

pyscf.pbc.tools.pbc.ifft(g, mesh)[source]

Perform the 3D inverse FFT from reciprocal (G) space to real (R) space.

Inverse FFT normalization factor is 1./N, same as in numpy.fft but different from MH (they use 1.).

Args:
g(nx*ny*nz,) ndarray

The function to be inverse FFT’d, flattened to a 1D array corresponding to the index order of span3.

mesh(3,) ndarray of ints (= nx,ny,nz)

The number G-vectors along each direction.

Returns:
(nx*ny*nz,) ndarray

The inverse FFT 1D array in same index order as Gv (natural order of numpy.fft).

pyscf.pbc.tools.pbc.ifftk(g, mesh, expikr)[source]

Perform the 3D inverse FFT of f(k+G) into a function which is (periodic*e^{ikr}).

fk(r) = (1/Ng) sum_G fk(k+G) e^{i(k+G)r} = (1/Ng) sum_G [fk(k+G)e^{iGr}] e^{ikr}

pyscf.pbc.tools.pbc.madelung(cell, kpts)[source]
pyscf.pbc.tools.pbc.mesh_to_cutoff(a, mesh)[source]

Convert #grid points to KE cutoff

pyscf.pbc.tools.pbc.precompute_exx(cell, kpts)[source]
pyscf.pbc.tools.pbc.round_to_cell0(r, tol=1e-06)[source]

Round scaled coordinates to reference unit cell

pyscf.pbc.tools.pbc.super_cell(cell, ncopy, wrap_around=False)[source]

Create an ncopy[0] x ncopy[1] x ncopy[2] supercell of the input cell Note this function differs from cell_plus_imgs() that cell_plus_imgs creates images in both +/- direction.

Args:

cell : instance of Cell

ncopy : (3,) array

wrap_aroundbool

Put the original cell centered on the super cell. It has the effects corresponding to the parameter wrap_around of cell.make_kpts.

Returns:

supcell : instance of Cell

pyscf.pbc.tools.print_funcs module

pyscf.pbc.tools.print_funcs.print_mo_energy_occ(mf, mo_energy, mo_occ, is_uhf)[source]
pyscf.pbc.tools.print_funcs.print_mo_energy_occ_kpts(mf, mo_energy_kpts, mo_occ_kpts, is_uhf)[source]

pyscf.pbc.tools.pyscf_ase module

pyscf.pbc.tools.pywannier90 module

pyscf.pbc.tools.tril module

pyscf.pbc.tools.tril.tril_index(ki, kj)[source]
pyscf.pbc.tools.tril.unpack_tril(in_array, nkpts, kp, kq, kr, ks)[source]

Module contents