pyscf.pbc.gto package

Submodules

pyscf.pbc.gto.cell module

pyscf.pbc.gto.cell.C(**kwargs)

This is a shortcut to build up Cell object.

Examples:

>>> from pyscf.pbc import gto
>>> cell = gto.M(a=numpy.eye(3)*4, atom='He 1 1 1', basis='6-31g')
class pyscf.pbc.gto.cell.Cell(**kwargs)

Bases: pyscf.gto.mole.Mole

A Cell object holds the basic information of a crystal.

Attributes:
a(3,3) ndarray

Lattice primitive vectors. Each row represents a lattice vector Reciprocal lattice vectors are given by b1,b2,b3 = 2 pi inv(a).T

mesh(3,) list of ints

The number G-vectors along each direction. The default value is estimated based on precision

pseudodict or str

To define pseudopotential.

precisionfloat

To control Ewald sums and lattice sums accuracy

rcutfloat

Cutoff radius (unit Bohr) in lattice summation. The default value is estimated based on the required precision.

ke_cutofffloat

If set, defines a spherical cutoff of planewaves, with .5 * G**2 < ke_cutoff The default value is estimated based on precision

dimensionint

Periodic dimensions. Default is 3

low_dim_ft_typestr

For semi-empirical periodic systems, whether to calculate integrals at the non-PBC dimension using the sampled mesh grids in infinity vacuum (inf_vacuum) or truncated Coulomb potential (analytic_2d_1). Unless explicitly specified, analytic_2d_1 is used for 2D system and inf_vacuum is assumed for 1D and 0D.

** Following attributes (for experts) are automatically generated. **

ew_eta, ew_cutfloat

The Ewald ‘eta’ and ‘cut’ parameters. See get_ewald_params()

(See other attributes in Mole)

Examples:

>>> mol = Mole(atom='H^2 0 0 0; H 0 0 1.1', basis='sto3g')
>>> cl = Cell()
>>> cl.build(a='3 0 0; 0 3 0; 0 0 3', atom='C 1 1 1', basis='sto3g')
>>> print(cl.atom_symbol(0))
C
property Gv
ao2mo(mo_coeffs, intor='int2e', erifile=None, dataname='eri_mo', **kwargs)

Integral transformation for arbitrary orbitals and arbitrary integrals. See more detalied documentation in func:ao2mo.kernel.

Args:
mo_coeffs (an np array or a list of arrays)A matrix of orbital

coefficients if it is a numpy ndarray, or four sets of orbital coefficients, corresponding to the four indices of (ij|kl).

Kwargs:
erifile (str or h5py File or h5py Group object)The file/object

to store the transformed integrals. If not given, the return value is an array (in memory) of the transformed integrals.

datanamestr

Note this argument is effective if erifile is given. The dataset name in the erifile (ref the hierarchy of HDF5 format http://www.hdfgroup.org/HDF5/doc1.6/UG/09_Groups.html). By assigning different dataname, the existed integral file can be reused. If the erifile contains the specified dataname, the old integrals will be replaced by the new one under the key dataname.

intor (str)integral name Name of the 2-electron integral. Ref

to getints_by_shell() for the complete list of available 2-electron integral names

Returns:

An array of transformed integrals if erifile is not given. Otherwise, return the file/fileobject if erifile is assigned.

Examples:

>>> import pyscf
>>> mol = pyscf.M(atom='O 0 0 0; H 0 1 0; H 0 0 1', basis='sto3g')
>>> mo1 = numpy.random.random((mol.nao_nr(), 10))
>>> mo2 = numpy.random.random((mol.nao_nr(), 8))
>>> eri1 = mol.ao2mo(mo1)
>>> print(eri1.shape)
(55, 55)
>>> eri1 = mol.ao2mo(mo1, compact=False)
>>> print(eri1.shape)
(100, 100)
>>> eri1 = mol.ao2mo(eri, (mo1,mo2,mo2,mo2))
>>> print(eri1.shape)
(80, 36)
>>> eri1 = mol.ao2mo(eri, (mo1,mo2,mo2,mo2), erifile='water.h5')
bas_rcut(bas_id, precision=1e-08)

Estimate the largest distance between the function and its image to reach the precision in overlap

precision ~ int g(r-0) g(r-Rcut)

build(dump_input=True, parse_arg=True, a=None, mesh=None, ke_cutoff=None, precision=None, nimgs=None, ew_eta=None, ew_cut=None, pseudo=None, basis=None, h=None, dimension=None, rcut=None, ecp=None, low_dim_ft_type=None, *args, **kwargs)

Setup Mole molecule and Cell and initialize some control parameters. Whenever you change the value of the attributes of Cell, you need call this function to refresh the internal data of Cell.

Kwargs:
a(3,3) ndarray

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

mesh(3,) ndarray of ints

The number of positive G-vectors along each direction.

ke_cutofffloat

If set, defines a spherical cutoff of planewaves, with .5 * G**2 < ke_cutoff The default value is estimated based on precision

precisionfloat

To control Ewald sums and lattice sums accuracy

nimgs(3,) ndarray of ints

Number of repeated images in lattice summation to produce periodicity. This value can be estimated based on the required precision. It’s recommended NOT making changes to this value.

rcutfloat

Cutoff radius (unit Bohr) in lattice summation to produce periodicity. The value can be estimated based on the required precision. It’s recommended NOT making changes to this value.

ew_eta, ew_cutfloat

Parameters eta and cut to converge Ewald summation. See get_ewald_params()

pseudodict or str

To define pseudopotential.

ecpdict or str

To define ECP type pseudopotential.

h(3,3) ndarray

a.T. Deprecated

dimensionint

Default is 3

low_dim_ft_typestr

For semi-empirical periodic systems, whether to calculate integrals at the non-PBC dimension using the sampled mesh grids in infinity vacuum (inf_vacuum) or truncated Coulomb potential (analytic_2d_1). Unless explicitly specified, analytic_2d_1 is used for 2D system and inf_vacuum is assumed for 1D and 0D.

copy()

Deepcopy of the given Mole object

property drop_exponent
dumps()

Serialize Cell object to a JSON formatted str.

energy_nuc(ew_eta=None, ew_cut=None)

Perform real (R) and reciprocal (G) space Ewald sum for the energy.

Formulation of Martin, App. F2.

Returns:
float

The Ewald energy consisting of overlap, self, and G-space sum.

See Also:

pyscf.pbc.gto.get_ewald_params

eval_ao(eval_name, coords, comp=None, kpts=None, kpt=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate PBC-AO function value on the given grids,

Args:

eval_name : str

Function

Expression

“GTOval_sph”

sum_T exp(ik*T) |AO>

“GTOval_ip_sph”

nabla sum_T exp(ik*T) |AO>

“GTOval_cart”

sum_T exp(ik*T) |AO>

“GTOval_ip_cart”

nabla sum_T exp(ik*T) |AO>

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

coords2D array, shape (N,3)

The coordinates of the grids.

Kwargs:
shls_slice2-element list

(shl_start, shl_end). If given, only part of AOs (shl_start <= shell_id < shl_end) are evaluated. By default, all shells defined in cell will be evaluated.

non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling dft.gen_grid.make_mask()

outndarray

If provided, results are written into this array.

Returns:

A list of 2D (or 3D) arrays to hold the AO values on grids. Each element of the list corresponds to a k-point and it has the shape (N,nao) Or shape (*,N,nao).

Examples:

>>> cell = pbc.gto.M(a=numpy.eye(3)*4, atom='He 1 1 1', basis='6-31g')
>>> coords = cell.get_uniform_grids([20,20,20])
>>> kpts = cell.make_kpts([3,3,3])
>>> ao_value = cell.pbc_eval_gto("GTOval_sph", coords, kpts)
>>> len(ao_value)
27
>>> ao_value[0].shape
(100, 2)
>>> ao_value = cell.pbc_eval_gto("GTOval_ig_sph", coords, kpts, comp=3)
>>> print(ao_value.shape)
>>> len(ao_value)
27
>>> ao_value[0].shape
(3, 100, 2)
eval_gto(eval_name, coords, comp=None, kpts=None, kpt=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate PBC-AO function value on the given grids,

Args:

eval_name : str

Function

Expression

“GTOval_sph”

sum_T exp(ik*T) |AO>

“GTOval_ip_sph”

nabla sum_T exp(ik*T) |AO>

“GTOval_cart”

sum_T exp(ik*T) |AO>

“GTOval_ip_cart”

nabla sum_T exp(ik*T) |AO>

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

coords2D array, shape (N,3)

The coordinates of the grids.

Kwargs:
shls_slice2-element list

(shl_start, shl_end). If given, only part of AOs (shl_start <= shell_id < shl_end) are evaluated. By default, all shells defined in cell will be evaluated.

non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling dft.gen_grid.make_mask()

outndarray

If provided, results are written into this array.

Returns:

A list of 2D (or 3D) arrays to hold the AO values on grids. Each element of the list corresponds to a k-point and it has the shape (N,nao) Or shape (*,N,nao).

Examples:

>>> cell = pbc.gto.M(a=numpy.eye(3)*4, atom='He 1 1 1', basis='6-31g')
>>> coords = cell.get_uniform_grids([20,20,20])
>>> kpts = cell.make_kpts([3,3,3])
>>> ao_value = cell.pbc_eval_gto("GTOval_sph", coords, kpts)
>>> len(ao_value)
27
>>> ao_value[0].shape
(100, 2)
>>> ao_value = cell.pbc_eval_gto("GTOval_ig_sph", coords, kpts, comp=3)
>>> print(ao_value.shape)
>>> len(ao_value)
27
>>> ao_value[0].shape
(3, 100, 2)
property ew_cut
property ew_eta
ewald(ew_eta=None, ew_cut=None)

Perform real (R) and reciprocal (G) space Ewald sum for the energy.

Formulation of Martin, App. F2.

Returns:
float

The Ewald energy consisting of overlap, self, and G-space sum.

See Also:

pyscf.pbc.gto.get_ewald_params

exp_to_discard = None
format_basis(basis_tab)

Convert the input Cell.basis to the internal data format:

{ atom: (l, kappa, ((-exp, c_1, c_2, ..), nprim, nctr, ptr-exps, ptr-contraction-coeff)), ... }
Args:
basis_tabdict

Similar to Cell.basis, it cannot be a str

Returns:

Formated basis

Examples:

>>> pbc.format_basis({'H':'gth-szv'})
{'H': [[0,
    (8.3744350009, -0.0283380461),
    (1.8058681460, -0.1333810052),
    (0.4852528328, -0.3995676063),
    (0.1658236932, -0.5531027541)]]}
format_pseudo(pseudo_tab)

Convert the input Cell.pseudo (dict) to the internal data format:

{ atom: ( (nelec_s, nele_p, nelec_d, ...),
         rloc, nexp, (cexp_1, cexp_2, ..., cexp_nexp),
         nproj_types,
         (r1, nproj1, ( (hproj1[1,1], hproj1[1,2], ..., hproj1[1,nproj1]),
                        (hproj1[2,1], hproj1[2,2], ..., hproj1[2,nproj1]),
                        ...
                        (hproj1[nproj1,1], hproj1[nproj1,2], ...        ) )),
         (r2, nproj2, ( (hproj2[1,1], hproj2[1,2], ..., hproj2[1,nproj1]),
         ... ) )
         )
 ... }
Args:
pseudo_tabdict

Similar to Cell.pseudo (a dict), it cannot be a str

Returns:

Formatted pseudo

Examples:

>>> pbc.format_pseudo({'H':'gth-blyp', 'He': 'gth-pade'})
{'H': [[1],
    0.2, 2, [-4.19596147, 0.73049821], 0],
 'He': [[2],
    0.2, 2, [-9.1120234, 1.69836797], 0]}
from_ase(ase_atom)

Update cell based on given ase atom object

Examples:

>>> from ase.lattice import bulk
>>> cell.from_ase(bulk('C', 'diamond', a=LATTICE_CONST))
gen_uniform_grids(mesh=None, **kwargs)

Generate a uniform real-space grid consistent w/ samp thm; see MH (3.19).

Args:

cell : instance of Cell

Returns:
coords(ngx*ngy*ngz, 3) ndarray

The real-space grid point coordinates.

get_Gv(mesh=None, **kwargs)

Calculate three-dimensional G-vectors for the cell; see MH (3.8).

Indices along each direction go as [0…N-1, -N…-1] to follow FFT convention.

Args:

cell : instance of Cell

Returns:
Gv(ngrids, 3) ndarray of floats

The array of G-vectors.

get_Gv_weights(mesh=None, **kwargs)

Calculate G-vectors and weights.

Returns:
Gv(ngris, 3) ndarray of floats

The array of G-vectors.

get_SI(Gv=None)

Calculate the structure factor (0D, 1D, 2D, 3D) for all atoms; see MH (3.34).

Args:

cell : instance of Cell

Gv(N,3) array

G vectors

Returns:
SI(natm, ngrids) ndarray, dtype=np.complex128

The structure factor for each atom at each G-vector.

get_abs_kpts(scaled_kpts)

Get absolute k-points (in 1/Bohr), given “scaled” k-points in fractions of lattice vectors.

Args:

scaled_kpts : (nkpts, 3) ndarray of floats

Returns:

abs_kpts : (nkpts, 3) ndarray of floats

get_bounding_sphere(rcut)

Finds all the lattice points within a sphere of radius rcut.

Defines a parallelipiped given by -N_x <= n_x <= N_x, with x in [1,3] See Martin p. 85

Args:
rcutnumber

real space cut-off for interaction

Returns:

cut : ndarray of 3 ints defining N_x

get_ewald_params(precision=1e-08, mesh=None)

Choose a reasonable value of Ewald ‘eta’ and ‘cut’ parameters. eta^2 is the exponent coefficient of the model Gaussian charge for nucleus at R: frac{eta^3}{pi^1.5} e^{-eta^2 (r-R)^2}

Choice is based on largest G vector and desired relative precision.

The relative error in the G-space sum is given by

precision ~ 4pi Gmax^2 e^{(-Gmax^2)/(4 eta^2)}

which determines eta. Then, real-space cutoff is determined by (exp. factors only)

precision ~ erfc(eta*rcut) / rcut ~ e^{(-eta**2 rcut*2)}

Returns:
ew_eta, ew_cutfloat

The Ewald ‘eta’ and ‘cut’ parameters.

get_kpts(nks, wrap_around=False, with_gamma_point=True, scaled_center=None)

Given number of kpoints along x,y,z , generate kpoints

Args:

nks : (3,) ndarray

Kwargs:
wrap_aroundbool

To ensure all kpts are in first Brillouin zone.

with_gamma_pointbool

Whether to shift Monkhorst-pack grid to include gamma-point.

scaled_center(3,) array

Shift all points in the Monkhorst-pack grid to be centered on scaled_center, given as the zeroth index of the returned kpts. Scaled meaning that the k-points are scaled to a grid from [-1,1] x [-1,1] x [-1,1]

Returns:

kpts in absolute value (unit 1/Bohr). Gamma point is placed at the first place in the k-points list

Examples:

>>> cell.make_kpts((4,4,4))
get_lattice_Ls(nimgs=None, rcut=None, dimension=None, discard=True)

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

get_nimgs(precision=None)

Choose number of basis function images in lattice sums to include for given precision in overlap, using

precision ~ int r^l e^{-alpha r^2} (r-rcut)^l e^{-alpha (r-rcut)^2} ~ (rcut^2/(2alpha))^l e^{alpha/2 rcut^2}

where alpha is the smallest exponent in the basis. Note that assumes an isolated exponent in the middle of the box, so it adds one additional lattice vector to be safe.

get_scaled_kpts(abs_kpts)

Get scaled k-points, given absolute k-points in 1/Bohr.

Args:

abs_kpts : (nkpts, 3) ndarray of floats

Returns:

scaled_kpts : (nkpts, 3) ndarray of floats

get_uniform_grids(mesh=None, **kwargs)

Generate a uniform real-space grid consistent w/ samp thm; see MH (3.19).

Args:

cell : instance of Cell

Returns:
coords(ngx*ngy*ngz, 3) ndarray

The real-space grid point coordinates.

property gs
property h
has_ecp()

Whether pseudo potential is used in the system.

kernel(dump_input=True, parse_arg=True, a=None, mesh=None, ke_cutoff=None, precision=None, nimgs=None, ew_eta=None, ew_cut=None, pseudo=None, basis=None, h=None, dimension=None, rcut=None, ecp=None, low_dim_ft_type=None, *args, **kwargs)

Setup Mole molecule and Cell and initialize some control parameters. Whenever you change the value of the attributes of Cell, you need call this function to refresh the internal data of Cell.

Kwargs:
a(3,3) ndarray

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

mesh(3,) ndarray of ints

The number of positive G-vectors along each direction.

ke_cutofffloat

If set, defines a spherical cutoff of planewaves, with .5 * G**2 < ke_cutoff The default value is estimated based on precision

precisionfloat

To control Ewald sums and lattice sums accuracy

nimgs(3,) ndarray of ints

Number of repeated images in lattice summation to produce periodicity. This value can be estimated based on the required precision. It’s recommended NOT making changes to this value.

rcutfloat

Cutoff radius (unit Bohr) in lattice summation to produce periodicity. The value can be estimated based on the required precision. It’s recommended NOT making changes to this value.

ew_eta, ew_cutfloat

Parameters eta and cut to converge Ewald summation. See get_ewald_params()

pseudodict or str

To define pseudopotential.

ecpdict or str

To define ECP type pseudopotential.

h(3,3) ndarray

a.T. Deprecated

dimensionint

Default is 3

low_dim_ft_typestr

For semi-empirical periodic systems, whether to calculate integrals at the non-PBC dimension using the sampled mesh grids in infinity vacuum (inf_vacuum) or truncated Coulomb potential (analytic_2d_1). Unless explicitly specified, analytic_2d_1 is used for 2D system and inf_vacuum is assumed for 1D and 0D.

lattice_vectors()

Convert the primitive lattice vectors.

Return 3x3 array in which each row represents one direction of the lattice vectors (unit in Bohr)

loads(molstr)

Deserialize a str containing a JSON document to a Cell object.

loads_(molstr)
make_ecp_env(_atm, _ecp, pre_env=[])

Generate the input arguments _ecpbas for ECP integrals

make_kpts(nks, wrap_around=False, with_gamma_point=True, scaled_center=None)

Given number of kpoints along x,y,z , generate kpoints

Args:

nks : (3,) ndarray

Kwargs:
wrap_aroundbool

To ensure all kpts are in first Brillouin zone.

with_gamma_pointbool

Whether to shift Monkhorst-pack grid to include gamma-point.

scaled_center(3,) array

Shift all points in the Monkhorst-pack grid to be centered on scaled_center, given as the zeroth index of the returned kpts. Scaled meaning that the k-points are scaled to a grid from [-1,1] x [-1,1] x [-1,1]

Returns:

kpts in absolute value (unit 1/Bohr). Gamma point is placed at the first place in the k-points list

Examples:

>>> cell.make_kpts((4,4,4))
property mesh
property nelec
property nimgs
pack()

Pack the input args of Cell to a dict, which can be serialized with pickle

pbc_eval_ao(eval_name, coords, comp=None, kpts=None, kpt=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate PBC-AO function value on the given grids,

Args:

eval_name : str

Function

Expression

“GTOval_sph”

sum_T exp(ik*T) |AO>

“GTOval_ip_sph”

nabla sum_T exp(ik*T) |AO>

“GTOval_cart”

sum_T exp(ik*T) |AO>

“GTOval_ip_cart”

nabla sum_T exp(ik*T) |AO>

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

coords2D array, shape (N,3)

The coordinates of the grids.

Kwargs:
shls_slice2-element list

(shl_start, shl_end). If given, only part of AOs (shl_start <= shell_id < shl_end) are evaluated. By default, all shells defined in cell will be evaluated.

non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling dft.gen_grid.make_mask()

outndarray

If provided, results are written into this array.

Returns:

A list of 2D (or 3D) arrays to hold the AO values on grids. Each element of the list corresponds to a k-point and it has the shape (N,nao) Or shape (*,N,nao).

Examples:

>>> cell = pbc.gto.M(a=numpy.eye(3)*4, atom='He 1 1 1', basis='6-31g')
>>> coords = cell.get_uniform_grids([20,20,20])
>>> kpts = cell.make_kpts([3,3,3])
>>> ao_value = cell.pbc_eval_gto("GTOval_sph", coords, kpts)
>>> len(ao_value)
27
>>> ao_value[0].shape
(100, 2)
>>> ao_value = cell.pbc_eval_gto("GTOval_ig_sph", coords, kpts, comp=3)
>>> print(ao_value.shape)
>>> len(ao_value)
27
>>> ao_value[0].shape
(3, 100, 2)
pbc_eval_gto(eval_name, coords, comp=None, kpts=None, kpt=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate PBC-AO function value on the given grids,

Args:

eval_name : str

Function

Expression

“GTOval_sph”

sum_T exp(ik*T) |AO>

“GTOval_ip_sph”

nabla sum_T exp(ik*T) |AO>

“GTOval_cart”

sum_T exp(ik*T) |AO>

“GTOval_ip_cart”

nabla sum_T exp(ik*T) |AO>

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

coords2D array, shape (N,3)

The coordinates of the grids.

Kwargs:
shls_slice2-element list

(shl_start, shl_end). If given, only part of AOs (shl_start <= shell_id < shl_end) are evaluated. By default, all shells defined in cell will be evaluated.

non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling dft.gen_grid.make_mask()

outndarray

If provided, results are written into this array.

Returns:

A list of 2D (or 3D) arrays to hold the AO values on grids. Each element of the list corresponds to a k-point and it has the shape (N,nao) Or shape (*,N,nao).

Examples:

>>> cell = pbc.gto.M(a=numpy.eye(3)*4, atom='He 1 1 1', basis='6-31g')
>>> coords = cell.get_uniform_grids([20,20,20])
>>> kpts = cell.make_kpts([3,3,3])
>>> ao_value = cell.pbc_eval_gto("GTOval_sph", coords, kpts)
>>> len(ao_value)
27
>>> ao_value[0].shape
(100, 2)
>>> ao_value = cell.pbc_eval_gto("GTOval_ig_sph", coords, kpts, comp=3)
>>> print(ao_value.shape)
>>> len(ao_value)
27
>>> ao_value[0].shape
(3, 100, 2)
pbc_intor(intor, comp=None, hermi=0, kpts=None, kpt=None, shls_slice=None, **kwargs)

One-electron integrals with PBC.

\[\sum_T \int \mu(r) * [intor] * \nu (r-T) dr\]

See also Mole.intor

precision = 1e-08
property rcut
reciprocal_vectors(norm_to=6.283185307179586)
\[\begin{split}\begin{align} \mathbf{b_1} &= 2\pi \frac{\mathbf{a_2} \times \mathbf{a_3}}{\mathbf{a_1} \cdot (\mathbf{a_2} \times \mathbf{a_3})} \\ \mathbf{b_2} &= 2\pi \frac{\mathbf{a_3} \times \mathbf{a_1}}{\mathbf{a_2} \cdot (\mathbf{a_3} \times \mathbf{a_1})} \\ \mathbf{b_3} &= 2\pi \frac{\mathbf{a_1} \times \mathbf{a_2}}{\mathbf{a_3} \cdot (\mathbf{a_1} \times \mathbf{a_2})} \end{align}\end{split}\]
to_mol()

Return a Mole object using the same atoms and basis functions as the Cell object.

tot_electrons(nkpts=1)

Total number of electrons

unpack(moldic)

Convert the packed dict to a Cell object, to generate the input arguments for Cell object.

unpack_(moldic)
property vol
pyscf.pbc.gto.cell.M(**kwargs)

This is a shortcut to build up Cell object.

Examples:

>>> from pyscf.pbc import gto
>>> cell = gto.M(a=numpy.eye(3)*4, atom='He 1 1 1', basis='6-31g')
pyscf.pbc.gto.cell.bas_rcut(cell, bas_id, precision=1e-08)

Estimate the largest distance between the function and its image to reach the precision in overlap

precision ~ int g(r-0) g(r-Rcut)

pyscf.pbc.gto.cell.classify_ecp_pseudo(cell, ecp, pp)
pyscf.pbc.gto.cell.conc_cell(cell1, cell2)

Concatenate two Cell objects.

pyscf.pbc.gto.cell.copy(cell)

Deepcopy of the given Cell object

pyscf.pbc.gto.cell.dumps(cell)

Serialize Cell object to a JSON formatted str.

pyscf.pbc.gto.cell.energy_nuc(cell, ew_eta=None, ew_cut=None)

Perform real (R) and reciprocal (G) space Ewald sum for the energy.

Formulation of Martin, App. F2.

Returns:
float

The Ewald energy consisting of overlap, self, and G-space sum.

See Also:

pyscf.pbc.gto.get_ewald_params

pyscf.pbc.gto.cell.error_for_ke_cutoff(cell, ke_cutoff)
pyscf.pbc.gto.cell.estimate_ke_cutoff(cell, precision=1e-08)

Energy cutoff estimation

pyscf.pbc.gto.cell.ewald(cell, ew_eta=None, ew_cut=None)

Perform real (R) and reciprocal (G) space Ewald sum for the energy.

Formulation of Martin, App. F2.

Returns:
float

The Ewald energy consisting of overlap, self, and G-space sum.

See Also:

pyscf.pbc.gto.get_ewald_params

pyscf.pbc.gto.cell.format_basis(basis_tab)

Convert the input Cell.basis to the internal data format:

{ atom: (l, kappa, ((-exp, c_1, c_2, ..), nprim, nctr, ptr-exps, ptr-contraction-coeff)), ... }
Args:
basis_tabdict

Similar to Cell.basis, it cannot be a str

Returns:

Formated basis

Examples:

>>> pbc.format_basis({'H':'gth-szv'})
{'H': [[0,
    (8.3744350009, -0.0283380461),
    (1.8058681460, -0.1333810052),
    (0.4852528328, -0.3995676063),
    (0.1658236932, -0.5531027541)]]}
pyscf.pbc.gto.cell.format_pseudo(pseudo_tab)

Convert the input Cell.pseudo (dict) to the internal data format:

{ atom: ( (nelec_s, nele_p, nelec_d, ...),
         rloc, nexp, (cexp_1, cexp_2, ..., cexp_nexp),
         nproj_types,
         (r1, nproj1, ( (hproj1[1,1], hproj1[1,2], ..., hproj1[1,nproj1]),
                        (hproj1[2,1], hproj1[2,2], ..., hproj1[2,nproj1]),
                        ...
                        (hproj1[nproj1,1], hproj1[nproj1,2], ...        ) )),
         (r2, nproj2, ( (hproj2[1,1], hproj2[1,2], ..., hproj2[1,nproj1]),
         ... ) )
         )
 ... }
Args:
pseudo_tabdict

Similar to Cell.pseudo (a dict), it cannot be a str

Returns:

Formatted pseudo

Examples:

>>> pbc.format_pseudo({'H':'gth-blyp', 'He': 'gth-pade'})
{'H': [[1],
    0.2, 2, [-4.19596147, 0.73049821], 0],
 'He': [[2],
    0.2, 2, [-9.1120234, 1.69836797], 0]}
pyscf.pbc.gto.cell.gen_uniform_grids(cell, mesh=None, **kwargs)

Generate a uniform real-space grid consistent w/ samp thm; see MH (3.19).

Args:

cell : instance of Cell

Returns:
coords(ngx*ngy*ngz, 3) ndarray

The real-space grid point coordinates.

pyscf.pbc.gto.cell.get_Gv(cell, mesh=None, **kwargs)

Calculate three-dimensional G-vectors for the cell; see MH (3.8).

Indices along each direction go as [0…N-1, -N…-1] to follow FFT convention.

Args:

cell : instance of Cell

Returns:
Gv(ngrids, 3) ndarray of floats

The array of G-vectors.

pyscf.pbc.gto.cell.get_Gv_weights(cell, mesh=None, **kwargs)

Calculate G-vectors and weights.

Returns:
Gv(ngris, 3) ndarray of floats

The array of G-vectors.

pyscf.pbc.gto.cell.get_SI(cell, Gv=None)

Calculate the structure factor (0D, 1D, 2D, 3D) for all atoms; see MH (3.34).

Args:

cell : instance of Cell

Gv(N,3) array

G vectors

Returns:
SI(natm, ngrids) ndarray, dtype=np.complex128

The structure factor for each atom at each G-vector.

pyscf.pbc.gto.cell.get_bounding_sphere(cell, rcut)

Finds all the lattice points within a sphere of radius rcut.

Defines a parallelipiped given by -N_x <= n_x <= N_x, with x in [1,3] See Martin p. 85

Args:
rcutnumber

real space cut-off for interaction

Returns:

cut : ndarray of 3 ints defining N_x

pyscf.pbc.gto.cell.get_ewald_params(cell, precision=1e-08, mesh=None)

Choose a reasonable value of Ewald ‘eta’ and ‘cut’ parameters. eta^2 is the exponent coefficient of the model Gaussian charge for nucleus at R: frac{eta^3}{pi^1.5} e^{-eta^2 (r-R)^2}

Choice is based on largest G vector and desired relative precision.

The relative error in the G-space sum is given by

precision ~ 4pi Gmax^2 e^{(-Gmax^2)/(4 eta^2)}

which determines eta. Then, real-space cutoff is determined by (exp. factors only)

precision ~ erfc(eta*rcut) / rcut ~ e^{(-eta**2 rcut*2)}

Returns:
ew_eta, ew_cutfloat

The Ewald ‘eta’ and ‘cut’ parameters.

pyscf.pbc.gto.cell.get_nimgs(cell, precision=None)

Choose number of basis function images in lattice sums to include for given precision in overlap, using

precision ~ int r^l e^{-alpha r^2} (r-rcut)^l e^{-alpha (r-rcut)^2} ~ (rcut^2/(2alpha))^l e^{alpha/2 rcut^2}

where alpha is the smallest exponent in the basis. Note that assumes an isolated exponent in the middle of the box, so it adds one additional lattice vector to be safe.

pyscf.pbc.gto.cell.get_uniform_grids(cell, mesh=None, **kwargs)

Generate a uniform real-space grid consistent w/ samp thm; see MH (3.19).

Args:

cell : instance of Cell

Returns:
coords(ngx*ngy*ngz, 3) ndarray

The real-space grid point coordinates.

pyscf.pbc.gto.cell.intor_cross(intor, cell1, cell2, comp=None, hermi=0, kpts=None, kpt=None, shls_slice=None, **kwargs)

1-electron integrals from two cells like

\[\langle \mu | intor | \nu \rangle, \mu \in cell1, \nu \in cell2\]
pyscf.pbc.gto.cell.loads(cellstr)

Deserialize a str containing a JSON document to a Cell object.

pyscf.pbc.gto.cell.make_kpts(cell, nks, wrap_around=False, with_gamma_point=True, scaled_center=None)

Given number of kpoints along x,y,z , generate kpoints

Args:

nks : (3,) ndarray

Kwargs:
wrap_aroundbool

To ensure all kpts are in first Brillouin zone.

with_gamma_pointbool

Whether to shift Monkhorst-pack grid to include gamma-point.

scaled_center(3,) array

Shift all points in the Monkhorst-pack grid to be centered on scaled_center, given as the zeroth index of the returned kpts. Scaled meaning that the k-points are scaled to a grid from [-1,1] x [-1,1] x [-1,1]

Returns:

kpts in absolute value (unit 1/Bohr). Gamma point is placed at the first place in the k-points list

Examples:

>>> cell.make_kpts((4,4,4))
pyscf.pbc.gto.cell.make_pseudo_env(cell, _atm, _pseudo, pre_env=[])
pyscf.pbc.gto.cell.pack(cell)

Pack the input args of Cell to a dict, which can be serialized with pickle

pyscf.pbc.gto.cell.tot_electrons(cell, nkpts=1)

Total number of electrons

pyscf.pbc.gto.cell.unpack(celldic)

Convert the packed dict to a Cell object, to generate the input arguments for Cell object.

pyscf.pbc.gto.ecp module

Short range part of ECP under PBC

pyscf.pbc.gto.ecp.ecp_int(cell, kpts=None)

pyscf.pbc.gto.eval_gto module

pyscf.pbc.gto.eval_gto.eval_gto(cell, eval_name, coords, comp=None, kpts=None, kpt=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate PBC-AO function value on the given grids,

Args:

eval_name : str

Function

Expression

“GTOval_sph”

sum_T exp(ik*T) |AO>

“GTOval_ip_sph”

nabla sum_T exp(ik*T) |AO>

“GTOval_cart”

sum_T exp(ik*T) |AO>

“GTOval_ip_cart”

nabla sum_T exp(ik*T) |AO>

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

coords2D array, shape (N,3)

The coordinates of the grids.

Kwargs:
shls_slice2-element list

(shl_start, shl_end). If given, only part of AOs (shl_start <= shell_id < shl_end) are evaluated. By default, all shells defined in cell will be evaluated.

non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling dft.gen_grid.make_mask()

outndarray

If provided, results are written into this array.

Returns:

A list of 2D (or 3D) arrays to hold the AO values on grids. Each element of the list corresponds to a k-point and it has the shape (N,nao) Or shape (*,N,nao).

Examples:

>>> cell = pbc.gto.M(a=numpy.eye(3)*4, atom='He 1 1 1', basis='6-31g')
>>> coords = cell.get_uniform_grids([20,20,20])
>>> kpts = cell.make_kpts([3,3,3])
>>> ao_value = cell.pbc_eval_gto("GTOval_sph", coords, kpts)
>>> len(ao_value)
27
>>> ao_value[0].shape
(100, 2)
>>> ao_value = cell.pbc_eval_gto("GTOval_ig_sph", coords, kpts, comp=3)
>>> print(ao_value.shape)
>>> len(ao_value)
27
>>> ao_value[0].shape
(3, 100, 2)
pyscf.pbc.gto.eval_gto.pbc_eval_gto(cell, eval_name, coords, comp=None, kpts=None, kpt=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate PBC-AO function value on the given grids,

Args:

eval_name : str

Function

Expression

“GTOval_sph”

sum_T exp(ik*T) |AO>

“GTOval_ip_sph”

nabla sum_T exp(ik*T) |AO>

“GTOval_cart”

sum_T exp(ik*T) |AO>

“GTOval_ip_cart”

nabla sum_T exp(ik*T) |AO>

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

coords2D array, shape (N,3)

The coordinates of the grids.

Kwargs:
shls_slice2-element list

(shl_start, shl_end). If given, only part of AOs (shl_start <= shell_id < shl_end) are evaluated. By default, all shells defined in cell will be evaluated.

non0tab2D bool array

mask array to indicate whether the AO values are zero. The mask array can be obtained by calling dft.gen_grid.make_mask()

outndarray

If provided, results are written into this array.

Returns:

A list of 2D (or 3D) arrays to hold the AO values on grids. Each element of the list corresponds to a k-point and it has the shape (N,nao) Or shape (*,N,nao).

Examples:

>>> cell = pbc.gto.M(a=numpy.eye(3)*4, atom='He 1 1 1', basis='6-31g')
>>> coords = cell.get_uniform_grids([20,20,20])
>>> kpts = cell.make_kpts([3,3,3])
>>> ao_value = cell.pbc_eval_gto("GTOval_sph", coords, kpts)
>>> len(ao_value)
27
>>> ao_value[0].shape
(100, 2)
>>> ao_value = cell.pbc_eval_gto("GTOval_ig_sph", coords, kpts, comp=3)
>>> print(ao_value.shape)
>>> len(ao_value)
27
>>> ao_value[0].shape
(3, 100, 2)

Module contents