pyscf.pbc.gto package#

Subpackages#

Submodules#

pyscf.pbc.gto.cell module#

pyscf.pbc.gto.cell.C(*args, **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)[source]#

Bases: MoleBase

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.

use_loose_rcutbool

If set to True, a loose rcut determined by shell radius is used, which is usually accurate enough for pure DFT calculations; otherwise, a tight rcut determined by overlap integral is used. Default value is False. Has no effect if rcut is set manually.

use_particle_mesh_ewaldbool

If set to True, use particle-mesh Ewald to compute the nuclear repulsion. Default value is False, meaning to use classical Ewald summation.

space_group_symmetrybool

Whether to consider space group symmetry. Default is False.

symmorphicbool

Whether the lattice is symmorphic. If set to True, even if the lattice is non-symmorphic, only symmorphic space group symmetry will be considered. Default is False, meaning the space group is determined by the lattice symmetry to be symmorphic or non-symmorphic.

lattice_symmetryNone or pbc.symm.Symmetry instance

The object containing the lattice symmetry information. Default is None.

(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#
bas_rcut(bas_id, precision=None)#

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=False, a=None, mesh=None, ke_cutoff=None, precision=None, nimgs=None, h=None, dimension=None, rcut=None, low_dim_ft_type=None, space_group_symmetry=None, symmorphic=None, use_loose_rcut=None, use_particle_mesh_ewald=None, fractional=None, *args, **kwargs)[source]#

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 cell lattice vectors. Each row represents a lattice vector.

fractionalbool

Whether the atom postions are specified in fractional coordinates. The default value is False, which means the coordinates are interpreted as Cartesian coordinate.

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.

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.

space_group_symmetrybool

Whether to consider space group symmetry. Default is False.

symmorphicbool

Whether the lattice is symmorphic. If set to True, even if the lattice is non-symmorphic, only symmorphic space group symmetry will be considered. Default is False, meaning the space group is determined by the lattice symmetry to be symmorphic or non-symmorphic.

build_lattice_symmetry(check_mesh_symmetry=True)[source]#

Build cell.lattice_symmetry object.

Kwargs:
check_mesh_symmetrybool

For nonsymmorphic symmetry groups, cell.mesh may have lower symmetry than the lattice. In this case, if check_mesh_symmetry is True, the lower symmetry group will be used. Otherwise, if check_mesh_symmetry is False, the mesh grid will be modified to satisfy the higher symmetry. Default value is True.

Note:

This function modifies the attributes of cell.

cutoff_to_mesh(ke_cutoff)[source]#

Convert KE cutoff to FFT-mesh

Args:
ke_cutofffloat

KE energy cutoff in a.u.

Returns:

mesh : (3,) array

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, cutoff=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()

cutofffloat

AO values smaller than cutoff will be set to zero. The default cutoff threshold is ~1e-22 (defined in gto/grid_ao_drv.h)

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([10,10,10])
>>> 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
(1000, 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, 1000, 2)
eval_gto(eval_name, coords, comp=None, kpts=None, kpt=None, shls_slice=None, non0tab=None, ao_loc=None, cutoff=None, out=None)[source]#

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()

cutofffloat

AO values smaller than cutoff will be set to zero. The default cutoff threshold is ~1e-22 (defined in gto/grid_ao_drv.h)

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([10,10,10])
>>> 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
(1000, 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, 1000, 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_atom(atoms, origin=0, axes=None, unit='angstrom')[source]#

Convert the input Mole.atom to the internal data format. Including, changing the nuclear charge to atom symbol, converting the coordinates to AU, rotate and shift the molecule. If the atom is a string, it takes “;” and “n” for the mark to separate atoms; “,” and arbitrary length of blank space to separate the individual terms for an atom. Blank line will be ignored.

Args:
atomslist or str

the same to Mole.atom

Kwargs:
originndarray

new axis origin.

axesndarray

(new_x, new_y, new_z), new coordinates

unitstr or number

If unit is one of strings (B, b, Bohr, bohr, AU, au), the coordinates of the input atoms are the atomic unit; If unit is one of strings (A, a, Angstrom, angstrom, Ang, ang), the coordinates are in the unit of angstrom; If a number is given, the number are considered as the Bohr value (in angstrom), which should be around 0.53. Set unit=1 if wishing to preserve the unit of the coordinates.

Returns:
“atoms” in the internal format. The internal format is
atom = [[atom1, (x, y, z)],
[atom2, (x, y, z)],
[atomN, (x, y, z)]]

Examples:

>>> gto.format_atom('9,0,0,0; h@1 0 0 1', origin=(1,1,1))
[['F', [-1.0, -1.0, -1.0]], ['H@1', [-1.0, -1.0, 0.0]]]
>>> gto.format_atom(['9,0,0,0', (1, (0, 0, 1))], origin=(1,1,1))
[['F', [-1.0, -1.0, -1.0]], ['H', [-1, -1, 0]]]
from_ase(ase_atom)[source]#

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, wrap_around=True)#

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, mesh=None, atmlst=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

atmlstlist of ints, optional

Indices of atoms for which the structure factors are computed.

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

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

get_abs_kpts(scaled_kpts)[source]#

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 parallelepiped 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_enuc(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

get_ewald_params(precision=None, 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, space_group_symmetry=False, time_reversal_symmetry=False, **kwargs)#

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]

space_group_symmetrybool

Whether to consider space group symmetry

time_reversal_symmetrybool

Whether to consider time reversal symmetry

Returns:

kpts in absolute value (unit 1/Bohr). Gamma point is placed at the first place in the k-points list; instance of KPoints if k-point symmetry is considered

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.

Kwargs:
discard:

Drop less important Ls based on AO values on grid

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_atom_coords(a=None)[source]#

Get scaled atomic coordinates.

get_scaled_kpts(abs_kpts, kpts_in_ibz=True)[source]#

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

Args:

abs_kpts : (nkpts, 3) ndarray of floats or KPoints object

kpts_in_ibzbool

If True, return k-points in IBZ; otherwise, return k-points in BZ. Default value is True. This has effects only if abs_kpts is a KPoints object

Returns:

scaled_kpts : (nkpts, 3) ndarray of floats

get_uniform_grids(mesh=None, wrap_around=True)#

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#
kernel(dump_input=True, parse_arg=False, a=None, mesh=None, ke_cutoff=None, precision=None, nimgs=None, h=None, dimension=None, rcut=None, low_dim_ft_type=None, space_group_symmetry=None, symmorphic=None, use_loose_rcut=None, use_particle_mesh_ewald=None, fractional=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 cell lattice vectors. Each row represents a lattice vector.

fractionalbool

Whether the atom postions are specified in fractional coordinates. The default value is False, which means the coordinates are interpreted as Cartesian coordinate.

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.

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.

space_group_symmetrybool

Whether to consider space group symmetry. Default is False.

symmorphicbool

Whether the lattice is symmorphic. If set to True, even if the lattice is non-symmorphic, only symmorphic space group symmetry will be considered. Default is False, meaning the space group is determined by the lattice symmetry to be symmorphic or non-symmorphic.

lattice_vectors()[source]#

Convert the primitive lattice vectors.

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

classmethod loads(molstr)[source]#

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

loads_(molstr)[source]#

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

make_kpts(nks, wrap_around=False, with_gamma_point=True, scaled_center=None, space_group_symmetry=False, time_reversal_symmetry=False, **kwargs)#

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]

space_group_symmetrybool

Whether to consider space group symmetry

time_reversal_symmetrybool

Whether to consider time reversal symmetry

Returns:

kpts in absolute value (unit 1/Bohr). Gamma point is placed at the first place in the k-points list; instance of KPoints if k-point symmetry is considered

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, cutoff=None, out=None, Ls=None, rcut=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()

cutofffloat

AO values smaller than cutoff will be set to zero. The default cutoff threshold is ~1e-22 (defined in gto/grid_ao_drv.h)

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([10,10,10])
>>> 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
(1000, 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, 1000, 2)
pbc_eval_gto(eval_name, coords, comp=None, kpts=None, kpt=None, shls_slice=None, non0tab=None, ao_loc=None, cutoff=None, out=None, Ls=None, rcut=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()

cutofffloat

AO values smaller than cutoff will be set to zero. The default cutoff threshold is ~1e-22 (defined in gto/grid_ao_drv.h)

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([10,10,10])
>>> 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
(1000, 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, 1000, 2)
pbc_intor(intor, comp=None, hermi=0, kpts=None, kpt=None, shls_slice=None, **kwargs)[source]#

One-electron integrals with PBC.

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

See also Mole.intor

precision = 1e-08#
property rcut#
rcut_by_shells(precision=None, rcut=0, return_pgf_radius=False)#

Compute shell and primitive gaussian function radii.

reciprocal_vectors(norm_to=6.283185307179586)[source]#
\[\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}\]
symmetrize_mesh(mesh=None)[source]#
to_mol()[source]#

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

tot_electrons(nkpts=1)#

Total number of electrons

classmethod unpack(moldic)[source]#

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

unpack_(moldic)[source]#

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

property vol#
pyscf.pbc.gto.cell.M(*args, **kwargs)[source]#

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=None)[source]#

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.conc_cell(cell1, cell2)[source]#

Concatenate two Cell objects.

pyscf.pbc.gto.cell.dumps(cell)[source]#

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, omega=None)[source]#

Error estimation based on nuclear attraction integrals

pyscf.pbc.gto.cell.estimate_ke_cutoff(cell, precision=None)[source]#

Energy cutoff estimation for nuclear attraction integrals

pyscf.pbc.gto.cell.estimate_rcut(cell, precision=None)[source]#

Lattice-sum cutoff for the entire system

pyscf.pbc.gto.cell.ewald(cell, ew_eta=None, ew_cut=None)[source]#

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.gen_uniform_grids(cell, mesh=None, wrap_around=True)#

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)[source]#

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)[source]#

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, mesh=None, atmlst=None)[source]#

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

atmlstlist of ints, optional

Indices of atoms for which the structure factors are computed.

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)[source]#

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

Defines a parallelepiped 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=None, mesh=None)[source]#

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)[source]#

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, wrap_around=True)[source]#

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)[source]#

1-electron integrals from two cells like

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

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, space_group_symmetry=False, time_reversal_symmetry=False, **kwargs)[source]#

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]

space_group_symmetrybool

Whether to consider space group symmetry

time_reversal_symmetrybool

Whether to consider time reversal symmetry

Returns:

kpts in absolute value (unit 1/Bohr). Gamma point is placed at the first place in the k-points list; instance of KPoints if k-point symmetry is considered

Examples:

>>> cell.make_kpts((4,4,4))
pyscf.pbc.gto.cell.pack(cell)[source]#

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

pyscf.pbc.gto.cell.pgf_rcut(l, alpha, coeff, precision=1e-08, rcut=0, max_cycle=10, eps=0.001)[source]#

Estimate the cutoff radii of primitive Gaussian functions based on their values in real space: c*rcut^(l+2)*exp(-alpha*rcut^2) ~ precision.

pyscf.pbc.gto.cell.rcut_by_shells(cell, precision=None, rcut=0, return_pgf_radius=False)[source]#

Compute shell and primitive gaussian function radii.

pyscf.pbc.gto.cell.tot_electrons(cell, nkpts=1)[source]#

Total number of electrons

pyscf.pbc.gto.cell.unpack(celldic)[source]#

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)[source]#

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, cutoff=None, out=None, Ls=None, rcut=None)[source]#

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()

cutofffloat

AO values smaller than cutoff will be set to zero. The default cutoff threshold is ~1e-22 (defined in gto/grid_ao_drv.h)

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([10,10,10])
>>> 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
(1000, 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, 1000, 2)
pyscf.pbc.gto.eval_gto.get_lattice_Ls(cell, nimgs=None, rcut=None, dimension=None, discard=True)[source]#

Get lattice-sum vectors for eval_gto

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, cutoff=None, out=None, Ls=None, rcut=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()

cutofffloat

AO values smaller than cutoff will be set to zero. The default cutoff threshold is ~1e-22 (defined in gto/grid_ao_drv.h)

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([10,10,10])
>>> 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
(1000, 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, 1000, 2)

pyscf.pbc.gto.ewald_methods module#

pyscf.pbc.gto.ewald_methods.bspline(u, ng, n=4, deriv=0)[source]#
pyscf.pbc.gto.ewald_methods.ewald_nuc_grad(cell, ew_eta=None, ew_cut=None)[source]#
pyscf.pbc.gto.ewald_methods.particle_mesh_ewald(cell, ew_eta=None, ew_cut=None, order=10)[source]#
pyscf.pbc.gto.ewald_methods.particle_mesh_ewald_nuc_grad(cell, ew_eta=None, ew_cut=None, order=10)[source]#

pyscf.pbc.gto.neighborlist module#

class pyscf.pbc.gto.neighborlist.NeighborListOpt(cell)[source]#

Bases: object

build(cell=None, cell1=None, Ls=None, ish_rcut=None, jsh_rcut=None, hermi=0, precision=None, set_nl=True, set_optimizer=True)[source]#
reset(free_nl=True)[source]#
pyscf.pbc.gto.neighborlist.build_neighbor_list_for_shlpairs(cell, cell1=None, Ls=None, ish_rcut=None, jsh_rcut=None, hermi=0, precision=None)[source]#

Build the neighbor list of shell pairs for periodic calculations.

Arguments:
cellpbc.gto.cell.Cell

The Cell instance for the bra basis functions.

cell1pbc.gto.cell.Cell, optional

The Cell instance for the ket basis functions. If not given, both bra and ket basis functions come from cell.

Ls(*,3) array, optional

The cartesian coordinates of the periodic images. Default is calculated by cell.get_lattice_Ls().

ish_rcut(nish,) array, optional

The cutoff radii of the shells for bra basis functions.

jsh_rcut(njsh,) array, optional

The cutoff radii of the shells for ket basis functions.

hermiint, optional

If \(hermi=1\), the task list is built only for the upper triangle of the matrix. Default is 0.

precisionfloat, optional

The integral precision. Default is cell.precision. If both ish_rcut and jsh_rcut are given, precision will be ignored.

Returns: ctypes.POINTER

The C pointer of the NeighborList structure.

pyscf.pbc.gto.neighborlist.free_neighbor_list(nl)[source]#
pyscf.pbc.gto.neighborlist.neighbor_list_to_ndarray(cell, cell1, nl)[source]#
Returns:
Ls_list: (nLtot,) ndarray

indices of Ls

Ls_idx: (2 x nish x njsh,) ndarray

starting and ending indices in Ls_list

Module contents#