pyscf.gto package

Submodules

pyscf.gto.cmd_args module

pyscf.gto.cmd_args.cmd_args()

get input from cmdline

pyscf.gto.ecp module

Effective core potential (ECP)

This module exposes some ecp integration functions from the C implementation.

Reference for ecp integral computation * Analytical integration J. Chem. Phys. 65, 3826 J. Chem. Phys. 111, 8778 J. Comput. Phys. 44, 289

  • Numerical integration

J. Comput. Chem. 27, 1009 Chem. Phys. Lett. 296, 445

pyscf.gto.ecp.core_configuration(nelec_core)
pyscf.gto.ecp.so_by_shell(mol, shls)

Spin-orbit coupling ECP in spinor basis i/2 <Pauli_matrix dot l U(r)>

pyscf.gto.ecp.type1_by_shell(mol, shls, cart=False)
pyscf.gto.ecp.type2_by_shell(mol, shls, cart=False)

pyscf.gto.eval_gto module

pyscf.gto.eval_gto.eval_gto(mol, eval_name, coords, comp=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate AO function value on the given grids,

Args:

eval_name : str

Function

comp

Expression

“GTOval_sph”

1

|AO>

“GTOval_ip_sph”

3

nabla |AO>

“GTOval_ig_sph”

3

(#C(0 1) g) |AO>

“GTOval_ipig_sph”

9

(#C(0 1) nabla g) |AO>

“GTOval_cart”

1

|AO>

“GTOval_ip_cart”

3

nabla |AO>

“GTOval_ig_cart”

3

(#C(0 1) g)|AO>

“GTOval_sph_deriv1”

4

GTO value and 1st order GTO values

“GTOval_sph_deriv2”

10

All derivatives up to 2nd order

“GTOval_sph_deriv3”

20

All derivatives up to 3rd order

“GTOval_sph_deriv4”

35

All derivatives up to 4th order

“GTOval_sp_spinor”

1

sigma dot p |AO> (spinor basis)

“GTOval_ipsp_spinor”

3

nabla sigma dot p |AO> (spinor basis)

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:
compint

Number of the components of the operator

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 mol 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:

2D array of shape (N,nao) Or 3D array of shape (*,N,nao) to store AO values on grids.

Examples:

>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz')
>>> coords = numpy.random.random((100,3))  # 100 random points
>>> ao_value = mol.eval_gto("GTOval_sph", coords)
>>> print(ao_value.shape)
(100, 24)
>>> ao_value = mol.eval_gto("GTOval_ig_sph", coords)
>>> print(ao_value.shape)
(3, 100, 24)

pyscf.gto.ft_ao module

Analytic Fourier transformation for AO and AO-pair value

pyscf.gto.ft_ao.ft_ao(mol, Gv, shls_slice=None, b=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), gxyz=None, Gvbase=None, verbose=None)

FT transform AO

pyscf.gto.ft_ao.ft_aopair(mol, Gv, shls_slice=None, aosym='s1', b=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), gxyz=None, Gvbase=None, buf=None, intor='GTO_ft_ovlp', comp=1, verbose=None)

FT transform AO pair int i(r) j(r) exp(-ikr) dr^3

pyscf.gto.mole module

Mole class and helper functions to handle paramters and attributes for GTO integrals. This module serves the interface to the integral library libcint.

pyscf.gto.mole.M(**kwargs)

This is a shortcut to build up Mole object.

Args: Same to Mole.build()

Examples:

>>> from pyscf import gto
>>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='6-31g')
class pyscf.gto.mole.Mole(**kwargs)

Bases: pyscf.lib.misc.StreamObject

Basic class to hold molecular structure and global options

Attributes:
verboseint

Print level

outputstr or None

Output file, default is None which dumps msg to sys.stdout

max_memoryint, float

Allowed memory in MB

chargeint

Charge of molecule. It affects the electron numbers

spinint or None

2S, num. alpha electrons - num. beta electrons to control multiplicity. If spin = None is set, multiplicity will be guessed based on the neutral molecule.

symmetrybool or str

Whether to use symmetry. When this variable is set to True, the molecule will be rotated and the highest rotation axis will be placed z-axis. If a string is given as the name of point group, the given point group symmetry will be used. Note that the input molecular coordinates will not be changed in this case.

symmetry_subgroupstr

subgroup

atomlist or str

To define molecluar structure. The internal format is

atom = [[atom1, (x, y, z)],
[atom2, (x, y, z)],
[atomN, (x, y, z)]]
unitstr

Angstrom or Bohr

basisdict or str

To define basis set.

nucmoddict or str or [function(nuc_charge, nucprop) => zeta]

Nuclear model. 0 or None means point nuclear model. Other values will enable Gaussian nuclear model. If a function is assigned to this attribute, the function will be called to generate the nuclear charge distribution value “zeta” and the relevant nuclear model will be set to Gaussian model. Default is point nuclear model.

nucpropdict

Nuclear properties (like g-factor ‘g’, quadrupole moments ‘Q’). It is needed by pyscf.prop module and submodules.

cartboolean

Using Cartesian GTO basis and integrals (6d,10f,15g)

** Following attributes are generated by Mole.build() **

stdoutfile object

Default is sys.stdout if Mole.output is not set

topgroupstr

Point group of the system.

groupnamestr

The supported subgroup of the point group. It can be one of Dooh, Coov, D2h, C2h, C2v, D2, Cs, Ci, C2, C1

nelectronint

sum of nuclear charges - Mole.charge

symm_orba list of numpy.ndarray

Symmetry adapted basis. Each element is a set of symm-adapted orbitals for one irreducible representation. The list index does not correspond to the id of irreducible representation.

irrep_ida list of int

Each element is one irreducible representation id associated with the basis stored in symm_orb. One irrep id stands for one irreducible representation symbol. The irrep symbol and the relevant id are defined in symm.param.IRREP_ID_TABLE

irrep_namea list of str

Each element is one irreducible representation symbol associated with the basis stored in symm_orb. The irrep symbols are defined in symm.param.IRREP_ID_TABLE

_builtbool

To label whether Mole.build() has been called. It is to ensure certain functions being initialized only once.

_basisdict

like Mole.basis, the internal format which is returned from the parser format_basis()

_keysa set of str

Store the keys appeared in the module. It is used to check misinput attributes

** Following attributes are arguments used by libcint library **

_atm :

[[charge, ptr-of-coord, nuc-model, ptr-zeta, 0, 0], [...]] each element reperesents one atom

natm :

number of atoms

_bas :

[[atom-id, angular-momentum, num-primitive-GTO, num-contracted-GTO, 0, ptr-of-exps, ptr-of-contract-coeff, 0], [...]] each element reperesents one shell

nbas :

number of shells

_env :

list of floats to store the coordinates, GTO exponents, contract-coefficients

Examples:

>>> mol = Mole(atom='H^2 0 0 0; H 0 0 1.1', basis='sto3g').build()
>>> print(mol.atom_symbol(0))
H^2
>>> print(mol.atom_pure_symbol(0))
H
>>> print(mol.nao_nr())
2
>>> print(mol.intor('int1e_ovlp_sph'))
[[ 0.99999999  0.43958641]
 [ 0.43958641  0.99999999]]
>>> mol.charge = 1
>>> mol.build()
<class 'pyscf.gto.mole.Mole'> has no attributes Charge
ao2mo(mo_coeffs, erifile=None, dataname='eri_mo', intor='int2e', **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')
ao_labels(fmt=True, base=0)

Labels of AO basis functions

Kwargs:
fmtstr or bool

if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-AO-notation)] or formatted strings based on the argument “fmt”

property ao_loc

Offset of every shell in the spherical basis function spectrum

Returns:

list, each entry is the corresponding start basis function id

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.ao_loc_nr(mol)
[0, 1, 2, 3, 6, 9, 10, 11, 12, 15, 18]
ao_loc_2c()

Offset of every shell in the spinor basis function spectrum

Returns:

list, each entry is the corresponding start id of spinor function

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.ao_loc_2c(mol)
[0, 2, 4, 6, 12, 18, 20, 22, 24, 30, 36]
ao_loc_nr(cart=None)

Offset of every shell in the spherical basis function spectrum

Returns:

list, each entry is the corresponding start basis function id

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.ao_loc_nr(mol)
[0, 1, 2, 3, 6, 9, 10, 11, 12, 15, 18]
aoslice_2c_by_atom()

2-component AO offset for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

aoslice_by_atom(ao_loc=None)

AO offsets for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

aoslice_nr_by_atom(ao_loc=None)

AO offsets for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

apply(fn, *args, **kwargs)

Apply the fn to rest arguments: return fn(*args, **kwargs). The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.

atom_charge(atm_id)

Nuclear effective charge of the given atom id Note “atom_charge /= charge(atom_symbol)” when ECP is enabled. Number of electrons screened by ECP can be obtained by charge(atom_symbol)-atom_charge

Args:
atm_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1')
>>> mol.atom_charge(1)
17
atom_charges()

np.asarray([mol.atom_charge(i) for i in range(mol.natm)])

atom_coord(atm_id, unit='Bohr')

Coordinates (ndarray) of the given atom id

Args:
atm_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1')
>>> mol.atom_coord(1)
[ 0.          0.          2.07869874]
atom_coords(unit='Bohr')

np.asarray([mol.atom_coords(i) for i in range(mol.natm)])

atom_mass_list(isotope_avg=False)

A list of mass for all atoms in the molecule

Kwargs:
isotope_avgboolean

Whether to use the isotope average mass as the atomic mass

atom_nelec_core(atm_id)

Number of core electrons for pseudo potential.

atom_nshells(atm_id)

Number of basis/shells of the given atom

Args:
atm_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1')
>>> mol.atom_nshells(1)
5
atom_pure_symbol(atm_id)

For the given atom id, return the standard symbol (striping special characters)

Args:
atm_idint

0-based

Examples:

>>> mol.build(atom='H^2 0 0 0; H 0 0 1.1')
>>> mol.atom_symbol(0)
H
atom_shell_ids(atm_id)

A list of the shell-ids of the given atom

Args:
atm_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1', basis='cc-pvdz')
>>> mol.atom_shell_ids(1)
[3, 4, 5, 6, 7]
atom_symbol(atm_id)

For the given atom id, return the input symbol (without striping special characters)

Args:
atm_idint

0-based

Examples:

>>> mol.build(atom='H^2 0 0 0; H 0 0 1.1')
>>> mol.atom_symbol(0)
H^2
bas_angular(bas_id)

The angular momentum associated with the given basis

Args:
bas_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1', basis='cc-pvdz')
>>> mol.bas_atom(7)
2
bas_atom(bas_id)

The atom (0-based id) that the given basis sits on

Args:
bas_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1', basis='cc-pvdz')
>>> mol.bas_atom(7)
1
bas_coord(bas_id)

Coordinates (ndarray) associated with the given basis id

Args:
bas_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1')
>>> mol.bas_coord(1)
[ 0.          0.          2.07869874]
bas_ctr_coeff(bas_id)

Contract coefficients (ndarray) of the given shell

Args:
bas_idint

0-based

Examples:

>>> mol.M(atom='H 0 0 0; Cl 0 0 1.1', basis='cc-pvdz')
>>> mol.bas_ctr_coeff(0)
[[ 10.03400444]
 [  4.1188704 ]
 [  1.53971186]]
bas_exp(bas_id)

exponents (ndarray) of the given shell

Args:
bas_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1', basis='cc-pvdz')
>>> mol.bas_kappa(0)
[ 13.01     1.962    0.4446]
bas_kappa(bas_id)

Kappa (if l < j, -l-1, else l) of the given shell

Args:
bas_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1', basis='cc-pvdz')
>>> mol.bas_kappa(3)
0
bas_len_cart(bas_id)

The number of Cartesian function associated with given basis

bas_len_spinor(bas_id)

The number of spinor associated with given basis If kappa is 0, return 4l+2

bas_nctr(bas_id)

The number of contracted GTOs for the given shell

Args:
bas_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1', basis='cc-pvdz')
>>> mol.bas_atom(3)
3
bas_nprim(bas_id)

The number of primitive GTOs for the given shell

Args:
bas_idint

0-based

Examples:

>>> mol.build(atom='H 0 0 0; Cl 0 0 1.1', basis='cc-pvdz')
>>> mol.bas_atom(3)
11
build(dump_input=True, parse_arg=True, verbose=None, output=None, max_memory=None, atom=None, basis=None, unit=None, nucmod=None, ecp=None, charge=None, spin=0, symmetry=None, symmetry_subgroup=None, cart=None)

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

Kwargs:
dump_inputbool

whether to dump the contents of input file in the output file

parse_argbool

whether to read the sys.argv and overwrite the relevant parameters

verboseint

Print level. If given, overwrite Mole.verbose

outputstr or None

Output file. If given, overwrite Mole.output

max_memoryint, float

Allowd memory in MB. If given, overwrite Mole.max_memory

atomlist or str

To define molecluar structure.

basisdict or str

To define basis set.

nucmoddict or str

Nuclear model. If given, overwrite Mole.nucmod

chargeint

Charge of molecule. It affects the electron numbers If given, overwrite Mole.charge

spinint

2S, num. alpha electrons - num. beta electrons to control multiplicity. If setting spin = None , multiplicity will be guessed based on the neutral molecule. If given, overwrite Mole.spin

symmetrybool or str

Whether to use symmetry. If given a string of point group name, the given point group symmetry will be used.

cart = False
cart2sph_coeff(normalized='sp')

Transformation matrix that transforms Cartesian GTOs to spherical GTOs for all basis functions

Kwargs:
normalizedstring or boolean

How the Cartesian GTOs are normalized. Except s and p functions, Cartesian GTOs do not have the universal normalization coefficients for the different components of the same shell. The value of this argument can be one of ‘sp’, ‘all’, None. ‘sp’ means the Cartesian s and p basis are normalized. ‘all’ means all Cartesian functions are normalized. None means none of the Cartesian functions are normalized. The default value ‘sp’ is the convention used by libcint library.

Examples:

>>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='ccpvtz')
>>> c = mol.cart2sph_coeff()
>>> s0 = mol.intor('int1e_ovlp_sph')
>>> s1 = c.T.dot(mol.intor('int1e_ovlp_cart')).dot(c)
>>> print(abs(s1-s0).sum())
>>> 4.58676826646e-15
cart_labels(fmt=True, base=0)

Labels of Cartesian GTO functions

Kwargs:

fmt : str or bool if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-xyz-notation)] or formatted strings based on the argument “fmt”

condense_to_shell(mat, compressor=<function amax>)

The given matrix is first partitioned to blocks, based on AO shell as delimiter. Then call compressor function to abstract each block.

copy()

Deepcopy of the given Mole object

dump_input()
dumps()

Serialize Mole object to a JSON formatted str.

property elements

A list of elements in the molecule

energy_nuc(charges=None, coords=None)

Compute nuclear repulsion energy (AU) or static Coulomb energy

Returns

float

etbs(etbs)

Generate even tempered basis. See also expand_etb()

Args:

etbs = [(l, n, alpha, beta), (l, n, alpha, beta),…]

Returns:

Formated basis

Examples:

>>> gto.expand_etbs([(0, 2, 1.5, 2.), (1, 2, 1, 2.)])
[[0, [6.0, 1]], [0, [3.0, 1]], [1, [1., 1]], [1, [2., 1]]]
eval_ao(eval_name, coords, comp=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate AO function value on the given grids,

Args:

eval_name : str

Function

comp

Expression

“GTOval_sph”

1

|AO>

“GTOval_ip_sph”

3

nabla |AO>

“GTOval_ig_sph”

3

(#C(0 1) g) |AO>

“GTOval_ipig_sph”

9

(#C(0 1) nabla g) |AO>

“GTOval_cart”

1

|AO>

“GTOval_ip_cart”

3

nabla |AO>

“GTOval_ig_cart”

3

(#C(0 1) g)|AO>

“GTOval_sph_deriv1”

4

GTO value and 1st order GTO values

“GTOval_sph_deriv2”

10

All derivatives up to 2nd order

“GTOval_sph_deriv3”

20

All derivatives up to 3rd order

“GTOval_sph_deriv4”

35

All derivatives up to 4th order

“GTOval_sp_spinor”

1

sigma dot p |AO> (spinor basis)

“GTOval_ipsp_spinor”

3

nabla sigma dot p |AO> (spinor basis)

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:
compint

Number of the components of the operator

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 mol 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:

2D array of shape (N,nao) Or 3D array of shape (*,N,nao) to store AO values on grids.

Examples:

>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz')
>>> coords = numpy.random.random((100,3))  # 100 random points
>>> ao_value = mol.eval_gto("GTOval_sph", coords)
>>> print(ao_value.shape)
(100, 24)
>>> ao_value = mol.eval_gto("GTOval_ig_sph", coords)
>>> print(ao_value.shape)
(3, 100, 24)
eval_gto(eval_name, coords, comp=None, shls_slice=None, non0tab=None, ao_loc=None, out=None)

Evaluate AO function value on the given grids,

Args:

eval_name : str

Function

comp

Expression

“GTOval_sph”

1

|AO>

“GTOval_ip_sph”

3

nabla |AO>

“GTOval_ig_sph”

3

(#C(0 1) g) |AO>

“GTOval_ipig_sph”

9

(#C(0 1) nabla g) |AO>

“GTOval_cart”

1

|AO>

“GTOval_ip_cart”

3

nabla |AO>

“GTOval_ig_cart”

3

(#C(0 1) g)|AO>

“GTOval_sph_deriv1”

4

GTO value and 1st order GTO values

“GTOval_sph_deriv2”

10

All derivatives up to 2nd order

“GTOval_sph_deriv3”

20

All derivatives up to 3rd order

“GTOval_sph_deriv4”

35

All derivatives up to 4th order

“GTOval_sp_spinor”

1

sigma dot p |AO> (spinor basis)

“GTOval_ipsp_spinor”

3

nabla sigma dot p |AO> (spinor basis)

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:
compint

Number of the components of the operator

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 mol 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:

2D array of shape (N,nao) Or 3D array of shape (*,N,nao) to store AO values on grids.

Examples:

>>> mol = gto.M(atom='O 0 0 0; H 0 0 1; H 0 1 0', basis='ccpvdz')
>>> coords = numpy.random.random((100,3))  # 100 random points
>>> ao_value = mol.eval_gto("GTOval_sph", coords)
>>> print(ao_value.shape)
(100, 24)
>>> ao_value = mol.eval_gto("GTOval_ig_sph", coords)
>>> print(ao_value.shape)
(3, 100, 24)
expand_etb(l, n, alpha, beta)

Generate the exponents of even tempered basis for Mole.basis. .. math:

e = e^{-\alpha * \beta^{i-1}} for i = 1 .. n
Args:
lint

Angular momentum

nint

Number of GTOs

Returns:

Formated basis

Examples:

>>> gto.expand_etb(1, 3, 1.5, 2)
[[1, [6.0, 1]], [1, [3.0, 1]], [1, [1.5, 1]]]
expand_etbs(etbs)

Generate even tempered basis. See also expand_etb()

Args:

etbs = [(l, n, alpha, beta), (l, n, alpha, beta),…]

Returns:

Formated basis

Examples:

>>> gto.expand_etbs([(0, 2, 1.5, 2.), (1, 2, 1, 2.)])
[[0, [6.0, 1]], [0, [3.0, 1]], [1, [1., 1]], [1, [2., 1]]]
format_atom(atom, origin=0, axes=None, unit='Ang')

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 spearate 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 reserve 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]]]
format_basis(basis_tab)

Convert the input Mole.basis to the internal data format.

``{ atom: [(l, ((-exp, c_1, c_2, ..),

(-exp, c_1, c_2, ..))),

(l, ((-exp, c_1, c_2, ..),

(-exp, c_1, c_2, ..)))], … }``

Args:
basis_tabdict

Similar to Mole.basis, it cannot be a str

Returns:

Formated basis

Examples:

>>> gto.format_basis({'H':'sto-3g', 'H^2': '3-21g'})
{'H': [[0,
    [3.4252509099999999, 0.15432897000000001],
    [0.62391373000000006, 0.53532813999999995],
    [0.16885539999999999, 0.44463454000000002]]],
 'H^2': [[0,
    [5.4471780000000001, 0.15628500000000001],
    [0.82454700000000003, 0.90469100000000002]],
    [0, [0.18319199999999999, 1.0]]]}
format_ecp(ecp_tab)

Convert the input ecp (dict) to the internal data format:

{ atom: (nelec,  # core electrons
((l, # l=-1 for UL, l>=0 for Ul to indicate |l><l|
(((exp_1, c_1), # for r^0

(exp_2, c_2), …),

((exp_1, c_1), # for r^1

(exp_2, c_2), …),

((exp_1, c_1), # for r^2

…))))),

…}

fromfile(filename, format=None)

Update the Mole object based on the input geometry file

fromstring(string, format='xyz')

Update the Mole object based on the input geometry string

get_enuc()
gto_norm(l, expnt)

Normalized factor for GTO radial part \(g=r^l e^{-\alpha r^2}\)

\[\frac{1}{\sqrt{\int g^2 r^2 dr}} = \sqrt{\frac{2^{2l+3} (l+1)! (2a)^{l+1.5}}{(2l+2)!\sqrt{\pi}}}\]

Ref: H. B. Schlegel and M. J. Frisch, Int. J. Quant. Chem., 54(1995), 83-87.

Args:
l (int):

angular momentum

expnt :

exponent \(\alpha\)

Returns:

normalization factor

Examples:

>>> print(gto_norm(0, 1))
2.5264751109842591
has_ecp()

Whether pseudo potential is used in the system.

incore_anyway = False
inertia_moment(mass=None, coords=None)
intor(intor, comp=None, hermi=0, aosym='s1', out=None, shls_slice=None)

Integral generator.

Args:
intorstr

Name of the 1e or 2e AO integrals. Ref to getints() for the complete list of available 1-electron integral names

Kwargs:
compint

Components of the integrals, e.g. int1e_ipovlp_sph has 3 components.

hermiint

Symmetry of the integrals

0 : no symmetry assumed (default)
1 : hermitian
2 : anti-hermitian
Returns:

ndarray of 1-electron integrals, can be either 2-dim or 3-dim, depending on comp

Examples:

>>> mol.build(atom='H 0 0 0; H 0 0 1.1', basis='sto-3g')
>>> mol.intor('int1e_ipnuc_sph', comp=3) # <nabla i | V_nuc | j>
[[[ 0.          0.        ]
  [ 0.          0.        ]]
 [[ 0.          0.        ]
  [ 0.          0.        ]]
 [[ 0.10289944  0.48176097]
  [-0.48176097 -0.10289944]]]
>>> mol.intor('int1e_nuc_spinor')
[[-1.69771092+0.j  0.00000000+0.j -0.67146312+0.j  0.00000000+0.j]
 [ 0.00000000+0.j -1.69771092+0.j  0.00000000+0.j -0.67146312+0.j]
 [-0.67146312+0.j  0.00000000+0.j -1.69771092+0.j  0.00000000+0.j]
 [ 0.00000000+0.j -0.67146312+0.j  0.00000000+0.j -1.69771092+0.j]]
intor_asymmetric(intor, comp=None)

One-electron integral generator. The integrals are assumed to be anti-hermitian

Args:
intorstr

Name of the 1-electron integral. Ref to getints() for the complete list of available 1-electron integral names

Kwargs:
compint

Components of the integrals, e.g. int1e_ipovlp has 3 components.

Returns:

ndarray of 1-electron integrals, can be either 2-dim or 3-dim, depending on comp

Examples:

>>> mol.build(atom='H 0 0 0; H 0 0 1.1', basis='sto-3g')
>>> mol.intor_asymmetric('int1e_nuc_spinor')
[[-1.69771092+0.j  0.00000000+0.j  0.67146312+0.j  0.00000000+0.j]
 [ 0.00000000+0.j -1.69771092+0.j  0.00000000+0.j  0.67146312+0.j]
 [-0.67146312+0.j  0.00000000+0.j -1.69771092+0.j  0.00000000+0.j]
 [ 0.00000000+0.j -0.67146312+0.j  0.00000000+0.j -1.69771092+0.j]]
intor_by_shell(intor, shells, comp=None)

For given 2, 3 or 4 shells, interface for libcint to get 1e, 2e, 2-center-2e or 3-center-2e integrals

Args:
intor_namestr

See also getints() for the supported intor_name

shlslist of int

The AO shell-ids of the integrals

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

Kwargs:
compint

Components of the integrals, e.g. int1e_ipovlp has 3 components.

Returns:

ndarray of 2-dim to 5-dim, depending on the integral type (1e, 2e, 3c-2e, 2c2e) and the value of comp

Examples:

The gradients of the spherical 2e integrals

>>> mol.build(atom='H 0 0 0; H 0 0 1.1', basis='sto-3g')
>>> gto.getints_by_shell('int2e_ip1_sph', (0,1,0,1), mol._atm, mol._bas, mol._env, comp=3)
[[[[[-0.        ]]]]
  [[[[-0.        ]]]]
  [[[[-0.08760462]]]]]
intor_symmetric(intor, comp=None)

One-electron integral generator. The integrals are assumed to be hermitian

Args:
intorstr

Name of the 1-electron integral. Ref to getints() for the complete list of available 1-electron integral names

Kwargs:
compint

Components of the integrals, e.g. int1e_ipovlp_sph has 3 components.

Returns:

ndarray of 1-electron integrals, can be either 2-dim or 3-dim, depending on comp

Examples:

>>> mol.build(atom='H 0 0 0; H 0 0 1.1', basis='sto-3g')
>>> mol.intor_symmetric('int1e_nuc_spinor')
[[-1.69771092+0.j  0.00000000+0.j -0.67146312+0.j  0.00000000+0.j]
 [ 0.00000000+0.j -1.69771092+0.j  0.00000000+0.j -0.67146312+0.j]
 [-0.67146312+0.j  0.00000000+0.j -1.69771092+0.j  0.00000000+0.j]
 [ 0.00000000+0.j -0.67146312+0.j  0.00000000+0.j -1.69771092+0.j]]
kernel(dump_input=True, parse_arg=True, verbose=None, output=None, max_memory=None, atom=None, basis=None, unit=None, nucmod=None, ecp=None, charge=None, spin=0, symmetry=None, symmetry_subgroup=None, cart=None)

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

Kwargs:
dump_inputbool

whether to dump the contents of input file in the output file

parse_argbool

whether to read the sys.argv and overwrite the relevant parameters

verboseint

Print level. If given, overwrite Mole.verbose

outputstr or None

Output file. If given, overwrite Mole.output

max_memoryint, float

Allowd memory in MB. If given, overwrite Mole.max_memory

atomlist or str

To define molecluar structure.

basisdict or str

To define basis set.

nucmoddict or str

Nuclear model. If given, overwrite Mole.nucmod

chargeint

Charge of molecule. It affects the electron numbers If given, overwrite Mole.charge

spinint

2S, num. alpha electrons - num. beta electrons to control multiplicity. If setting spin = None , multiplicity will be guessed based on the neutral molecule. If given, overwrite Mole.spin

symmetrybool or str

Whether to use symmetry. If given a string of point group name, the given point group symmetry will be used.

loads(molstr)

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

loads_(molstr)
make_atm_env(atom, ptr=0, nucmod=1, nucprop=None)

Convert the internal format Mole._atom to the format required by libcint integrals

make_bas_env(basis_add, atom_id=0, ptr=0)

Convert Mole.basis to the argument bas for libcint integrals

make_ecp_env(_atm, _ecp, pre_env=[])

Generate the input arguments _ecpbas for ECP integrals

make_env(atoms, basis, pre_env=[], nucmod={}, nucprop=None)

Generate the input arguments for libcint library based on internal format Mole._atom and Mole._basis

property ms

Spin quantum number. multiplicity = ms*2+1

property multiplicity
property nao
nao_2c()

Total number of contracted spinor GTOs for the given Mole object

nao_2c_range(bas_id0, bas_id1)

Lower and upper boundary of contracted spinor basis functions associated with the given shell range

Args:
mol :

Mole object

bas_id0int

start shell id, 0-based

bas_id1int

stop shell id, 0-based

Returns:

tupel of start basis function id and the stop function id

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.nao_2c_range(mol, 2, 4)
(4, 12)
nao_cart()

Total number of contracted cartesian GTOs for the given Mole object

nao_nr(cart=None)

Total number of contracted GTOs for the given Mole object

nao_nr_range(bas_id0, bas_id1)

Lower and upper boundary of contracted spherical basis functions associated with the given shell range

Args:
mol :

Mole object

bas_id0int

start shell id

bas_id1int

stop shell id

Returns:

tupel of start basis function id and the stop function id

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.nao_nr_range(mol, 2, 4)
(2, 6)
property natm
property nbas
property nelec
property nelectron
npgto_nr(cart=None)

Total number of primitive spherical GTOs for the given Mole object

offset_2c_by_atom()

2-component AO offset for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

offset_ao_by_atom(ao_loc=None)

AO offsets for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

offset_nr_by_atom(ao_loc=None)

AO offsets for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

property omega
pack()

Pack the input args of Mole to a dict.

Note this function only pack the input arguments (not the entire Mole class). Modifications to mol._atm, mol._bas, mol._env are not tracked. Use dumps() to serialize the entire Mole object.

search_ao_label(label)

Find the index of the AO basis function based on the given ao_label

Args:
ao_labelstring or a list of strings

The regular expression pattern to match the orbital labels returned by mol.ao_labels()

Returns:

A list of index for the AOs that matches the given ao_label RE pattern

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='ccpvtz')
>>> mol.search_ao_label('Cl.*p')
[19 20 21 22 23 24 25 26 27 28 29 30]
>>> mol.search_ao_label('Cl 2p')
[19 20 21]
>>> mol.search_ao_label(['Cl.*d', 'Cl 4p'])
[25 26 27 31 32 33 34 35 36 37 38 39 40]
search_ao_nr(atm_id, l, m, atmshell)

Search the first basis function id (not the shell id) which matches the given atom-id, angular momentum magnetic angular momentum, principal shell.

Args:
atm_idint

atom id, 0-based

lint

angular momentum

mint

magnetic angular momentum

atmshellint

principal quantum number

Returns:

basis function id, 0-based. If not found, return None

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='sto-3g')
>>> mol.search_ao_nr(1, 1, -1, 3) # Cl 3px
7
search_ao_r(atm_id, l, j, m, atmshell)
search_shell_id(atm_id, l)
set_common_orig(coord)

Update common origin for integrals of dipole, rxp etc. Note the unit of the coordinates needs to be Bohr

Examples:

>>> mol.set_common_origin(0)
>>> mol.set_common_origin((1,0,0))
set_common_orig_(coord)

Update common origin for integrals of dipole, rxp etc. Note the unit of the coordinates needs to be Bohr

Examples:

>>> mol.set_common_origin(0)
>>> mol.set_common_origin((1,0,0))
set_common_origin(coord)

Update common origin for integrals of dipole, rxp etc. Note the unit of the coordinates needs to be Bohr

Examples:

>>> mol.set_common_origin(0)
>>> mol.set_common_origin((1,0,0))
set_common_origin_(coord)

Update common origin for integrals of dipole, rxp etc. Note the unit of the coordinates needs to be Bohr

Examples:

>>> mol.set_common_origin(0)
>>> mol.set_common_origin((1,0,0))
set_f12_zeta(zeta)

Set zeta for YP exp(-zeta r12)/r12 or STG exp(-zeta r12) type integrals

set_geom_(atoms_or_coords, unit=None, symmetry=None, inplace=True)

Update geometry

set_nuc_mod(atm_id, zeta)

Change the nuclear charge distribution of the given atom ID. The charge distribution is defined as: rho(r) = nuc_charge * Norm * exp(-zeta * r^2). This function can only be called after .build() method is executed.

Examples:

>>> for ia in range(mol.natm):
...     zeta = gto.filatov_nuc_mod(mol.atom_charge(ia))
...     mol.set_nuc_mod(ia, zeta)
set_nuc_mod_(atm_id, zeta)

Change the nuclear charge distribution of the given atom ID. The charge distribution is defined as: rho(r) = nuc_charge * Norm * exp(-zeta * r^2). This function can only be called after .build() method is executed.

Examples:

>>> for ia in range(mol.natm):
...     zeta = gto.filatov_nuc_mod(mol.atom_charge(ia))
...     mol.set_nuc_mod(ia, zeta)
set_range_coulomb(omega)

Switch on range-separated Coulomb operator for all 2e integrals

Args:

omega : double

= 0 : Regular electron repulsion integral
> 0 : Long-range operator erf(omega r12) / r12
< 0 : Short-range operator erfc(omega r12) /r12
set_range_coulomb_(omega)

Switch on range-separated Coulomb operator for all 2e integrals

Args:

omega : double

= 0 : Regular electron repulsion integral
> 0 : Long-range operator erf(omega r12) / r12
< 0 : Short-range operator erfc(omega r12) /r12
set_rinv_orig(coord)

Update origin for operator \(\frac{1}{|r-R_O|}\). Note the unit is Bohr

Examples:

>>> mol.set_rinv_origin(0)
>>> mol.set_rinv_origin((0,1,0))
set_rinv_orig_(coord)

Update origin for operator \(\frac{1}{|r-R_O|}\). Note the unit is Bohr

Examples:

>>> mol.set_rinv_origin(0)
>>> mol.set_rinv_origin((0,1,0))
set_rinv_origin(coord)

Update origin for operator \(\frac{1}{|r-R_O|}\). Note the unit is Bohr

Examples:

>>> mol.set_rinv_origin(0)
>>> mol.set_rinv_origin((0,1,0))
set_rinv_origin_(coord)

Update origin for operator \(\frac{1}{|r-R_O|}\). Note the unit is Bohr

Examples:

>>> mol.set_rinv_origin(0)
>>> mol.set_rinv_origin((0,1,0))
set_rinv_zeta(zeta)

Assume the charge distribution on the “rinv_origin”. zeta is the parameter to control the charge distribution: rho(r) = Norm * exp(-zeta * r^2). Be careful when call this function. It affects the behavior of int1e_rinv_* functions. Make sure to set it back to 0 after using it!

set_rinv_zeta_(zeta)

Assume the charge distribution on the “rinv_origin”. zeta is the parameter to control the charge distribution: rho(r) = Norm * exp(-zeta * r^2). Be careful when call this function. It affects the behavior of int1e_rinv_* functions. Make sure to set it back to 0 after using it!

sph2spinor_coeff()

Transformation matrix that transforms real-spherical GTOs to spinor GTOs for all basis functions

Examples:

>>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='ccpvtz')
>>> ca, cb = mol.sph2spinor_coeff()
>>> s0 = mol.intor('int1e_ovlp_spinor')
>>> s1 = ca.conj().T.dot(mol.intor('int1e_ovlp_sph')).dot(ca)
>>> s1+= cb.conj().T.dot(mol.intor('int1e_ovlp_sph')).dot(cb)
>>> print(abs(s1-s0).max())
>>> 6.66133814775e-16
sph_labels(fmt=True, base=0)

Labels for spherical GTO functions

Kwargs:

fmt : str or bool if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-real-spherical-notation] or formatted strings based on the argument “fmt”

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='sto-3g')
>>> gto.sph_labels(mol)
[(0, 'H', '1s', ''), (1, 'Cl', '1s', ''), (1, 'Cl', '2s', ''), (1, 'Cl', '3s', ''),
 (1, 'Cl', '2p', 'x'), (1, 'Cl', '2p', 'y'), (1, 'Cl', '2p', 'z'), (1, 'Cl', '3p', 'x'),
 (1, 'Cl', '3p', 'y'), (1, 'Cl', '3p', 'z')]
spheric_labels(fmt=True, base=0)

Labels for spherical GTO functions

Kwargs:

fmt : str or bool if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-real-spherical-notation] or formatted strings based on the argument “fmt”

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='sto-3g')
>>> gto.sph_labels(mol)
[(0, 'H', '1s', ''), (1, 'Cl', '1s', ''), (1, 'Cl', '2s', ''), (1, 'Cl', '3s', ''),
 (1, 'Cl', '2p', 'x'), (1, 'Cl', '2p', 'y'), (1, 'Cl', '2p', 'z'), (1, 'Cl', '3p', 'x'),
 (1, 'Cl', '3p', 'y'), (1, 'Cl', '3p', 'z')]
spinor_labels(fmt=True, base=0)

Labels of spinor GTO functions

time_reversal_map()

The index to map the spinor functions and its time reversal counterpart. The returned indices have postive or negative values. For the i-th basis function, if the returned j = idx[i] < 0, it means \(T|i\rangle = -|j\rangle\), otherwise \(T|i\rangle = |j\rangle\)

tmap()

The index to map the spinor functions and its time reversal counterpart. The returned indices have postive or negative values. For the i-th basis function, if the returned j = idx[i] < 0, it means \(T|i\rangle = -|j\rangle\), otherwise \(T|i\rangle = |j\rangle\)

to_uncontracted_cartesian_basis()

Decontract the basis of a Mole or a Cell. Returns a Mole (Cell) object with the uncontracted basis environment and a list of coefficients that transform the uncontracted cartesian basis to the original basis. Each element in the list corresponds to one shell of the original Mole (Cell).

Examples:

>>> mol = gto.M(atom='Ne', basis='ccpvdz')
>>> pmol, ctr_coeff = mol.to_uncontracted_cartesian_basis()
>>> c = scipy.linalg.block_diag(*ctr_coeff)
>>> s = reduce(numpy.dot, (c.T, pmol.intor('int1e_ovlp'), c))
>>> abs(s-mol.intor('int1e_ovlp')).max()
0.0
tofile(filename, format=None)

Write molecular geometry to a file of the required format.

Supported output formats:
raw: Each line is <symobl> <x> <y> <z>
xyz: XYZ cartesian coordinates format
zmat: Z-matrix format
tostring(format='raw')

Convert molecular geometry to a string of the required format.

Supported output formats:
raw: Each line is <symobl> <x> <y> <z>
xyz: XYZ cartesian coordinates format
zmat: Z-matrix format
tot_electrons()

Total number of electrons for the given molecule

Returns:

electron number in integer

Examples:

>>> mol = gto.M(atom='H 0 1 0; C 0 0 1', charge=1)
>>> mol.tot_electrons()
6
unit = 'angstrom'
unpack(moldic)

Unpack a dict which is packed by pack(), to generate the input arguments for Mole object.

unpack_(moldic)
update(chkfile)
update_from_chk(chkfile)
verbose = 3
with_common_orig(coord)

Retuen a temporary mol context which has the rquired common origin. The required common origin has no effects out of the temporary context. See also mol.set_common_origin()

Examples:

>>> with mol.with_common_origin((1,0,0)):
...     mol.intor('int1e_r', comp=3)
with_common_origin(coord)

Retuen a temporary mol context which has the rquired common origin. The required common origin has no effects out of the temporary context. See also mol.set_common_origin()

Examples:

>>> with mol.with_common_origin((1,0,0)):
...     mol.intor('int1e_r', comp=3)
with_integral_screen(threshold)

Retuen a temporary mol context which has the rquired integral screen threshold

with_long_range_coulomb(omega)

Retuen a temporary mol context for long-range part of range-separated Coulomb operator.

with_range_coulomb(omega)

Retuen a temporary mol context which sets the required parameter omega for range-separated Coulomb operator. If omega = None, return the context for regular Coulomb integrals. See also mol.set_range_coulomb()

Examples:

>>> with mol.with_range_coulomb(omega=1.5):
...     mol.intor('int2e')
with_rinv_as_nucleus(atm_id)

Retuen a temporary mol context in which the rinv operator (1/r) is treated like the Coulomb potential of a Gaussian charge distribution rho(r) = Norm * exp(-zeta * r^2) at the place of the input atm_id.

Examples:

>>> with mol.with_rinv_at_nucleus(3):
...     mol.intor('int1e_rinv')
with_rinv_at_nucleus(atm_id)

Retuen a temporary mol context in which the rinv operator (1/r) is treated like the Coulomb potential of a Gaussian charge distribution rho(r) = Norm * exp(-zeta * r^2) at the place of the input atm_id.

Examples:

>>> with mol.with_rinv_at_nucleus(3):
...     mol.intor('int1e_rinv')
with_rinv_orig(coord)

Retuen a temporary mol context which has the rquired origin of 1/r operator. The required origin has no effects out of the temporary context. See also mol.set_rinv_origin()

Examples:

>>> with mol.with_rinv_origin((1,0,0)):
...     mol.intor('int1e_rinv')
with_rinv_origin(coord)

Retuen a temporary mol context which has the rquired origin of 1/r operator. The required origin has no effects out of the temporary context. See also mol.set_rinv_origin()

Examples:

>>> with mol.with_rinv_origin((1,0,0)):
...     mol.intor('int1e_rinv')
with_rinv_zeta(zeta)

Retuen a temporary mol context which has the rquired Gaussian charge distribution placed at “rinv_origin”: rho(r) = Norm * exp(-zeta * r^2). See also mol.set_rinv_zeta()

Examples:

>>> with mol.with_rinv_zeta(zeta=1.5), mol.with_rinv_origin((1.,0,0)):
...     mol.intor('int1e_rinv')
with_short_range_coulomb(omega)

Retuen a temporary mol context for short-range part of range-separated Coulomb operator.

pyscf.gto.mole.ao_labels(mol, fmt=True, base=0)

Labels of AO basis functions

Kwargs:
fmtstr or bool

if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-AO-notation)] or formatted strings based on the argument “fmt”

pyscf.gto.mole.ao_loc_2c(mol)

Offset of every shell in the spinor basis function spectrum

Returns:

list, each entry is the corresponding start id of spinor function

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.ao_loc_2c(mol)
[0, 2, 4, 6, 12, 18, 20, 22, 24, 30, 36]
pyscf.gto.mole.ao_loc_nr(mol, cart=None)

Offset of every shell in the spherical basis function spectrum

Returns:

list, each entry is the corresponding start basis function id

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.ao_loc_nr(mol)
[0, 1, 2, 3, 6, 9, 10, 11, 12, 15, 18]
pyscf.gto.mole.aoslice_by_atom(mol, ao_loc=None)

AO offsets for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

pyscf.gto.mole.atom_mass_list(mol, isotope_avg=False)

A list of mass for all atoms in the molecule

Kwargs:
isotope_avgboolean

Whether to use the isotope average mass as the atomic mass

pyscf.gto.mole.atom_types(atoms, basis=None)

symmetry inequivalent atoms

pyscf.gto.mole.cart2j_kappa(kappa, l=None, normalized=None)

Cartesian to spinor transformation matrix for kappa

Kwargs:
normalized :

How the Cartesian GTOs are normalized. ‘sp’ means the s and p functions are normalized (this is the convention used by libcint library).

pyscf.gto.mole.cart2j_l(l, normalized=None)

Cartesian to spinor transformation matrix for angular moment l

Kwargs:
normalized :

How the Cartesian GTOs are normalized. ‘sp’ means the s and p functions are normalized (this is the convention used by libcint library).

pyscf.gto.mole.cart2sph(l, c_tensor=None, normalized=None)

Cartesian to real spherical transformation matrix

Kwargs:
normalized :

How the Cartesian GTOs are normalized. ‘sp’ means the s and p functions are normalized (this is the convention used by libcint library).

pyscf.gto.mole.cart2spinor_kappa(kappa, l=None, normalized=None)

Cartesian to spinor transformation matrix for kappa

Kwargs:
normalized :

How the Cartesian GTOs are normalized. ‘sp’ means the s and p functions are normalized (this is the convention used by libcint library).

pyscf.gto.mole.cart2spinor_l(l, normalized=None)

Cartesian to spinor transformation matrix for angular moment l

Kwargs:
normalized :

How the Cartesian GTOs are normalized. ‘sp’ means the s and p functions are normalized (this is the convention used by libcint library).

pyscf.gto.mole.cart2zmat(coord)
>>> c = numpy.array((
(0.000000000000,  1.889726124565,  0.000000000000),
(0.000000000000,  0.000000000000, -1.889726124565),
(1.889726124565, -1.889726124565,  0.000000000000),
(1.889726124565,  0.000000000000,  1.133835674739)))
>>> print(cart2zmat(c))
1
1 2.67247631453057
1 4.22555607338457 2 50.7684795164077
1 2.90305235726773 2 79.3904651036893 3 6.20854462618583
pyscf.gto.mole.cart_labels(mol, fmt=True, base=0)

Labels of Cartesian GTO functions

Kwargs:

fmt : str or bool if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-xyz-notation)] or formatted strings based on the argument “fmt”

pyscf.gto.mole.chiral_mol(mol1, mol2=None)

Detect whether the given molelcule is chiral molecule or two molecules are chiral isomers.

pyscf.gto.mole.conc_env(atm1, bas1, env1, atm2, bas2, env2)

Concatenate two Mole’s integral parameters. This function can be used to construct the environment for cross integrals like

\[\langle \mu | \nu \rangle, \mu \in mol1, \nu \in mol2\]
Returns:

Concatenated atm, bas, env

Examples:

Compute the overlap between H2 molecule and O atom

>>> mol1 = gto.M(atom='H 0 1 0; H 0 0 1', basis='sto3g')
>>> mol2 = gto.M(atom='O 0 0 0', basis='sto3g')
>>> atm3, bas3, env3 = gto.conc_env(mol1._atm, mol1._bas, mol1._env,
...                                 mol2._atm, mol2._bas, mol2._env)
>>> gto.moleintor.getints('int1e_ovlp_sph', atm3, bas3, env3, range(2), range(2,5))
[[ 0.04875181  0.44714688  0.          0.37820346  0.        ]
 [ 0.04875181  0.44714688  0.          0.          0.37820346]]
pyscf.gto.mole.conc_mol(mol1, mol2)

Concatenate two Mole objects.

pyscf.gto.mole.condense_to_shell(mol, mat, compressor=<function amax>)

The given matrix is first partitioned to blocks, based on AO shell as delimiter. Then call compressor function to abstract each block.

pyscf.gto.mole.copy(mol)

Deepcopy of the given Mole object

pyscf.gto.mole.dumps(mol)

Serialize Mole object to a JSON formatted str.

pyscf.gto.mole.dyall_nuc_mod(nuc_charge, nucprop={})

Generate the nuclear charge distribution parameter zeta rho(r) = nuc_charge * Norm * exp(-zeta * r^2)

Ref. L. Visscher and K. Dyall, At. Data Nucl. Data Tables, 67, 207 (1997)

pyscf.gto.mole.energy_nuc(mol, charges=None, coords=None)

Compute nuclear repulsion energy (AU) or static Coulomb energy

Returns

float

pyscf.gto.mole.etbs(etbs)

Generate even tempered basis. See also expand_etb()

Args:

etbs = [(l, n, alpha, beta), (l, n, alpha, beta),…]

Returns:

Formated basis

Examples:

>>> gto.expand_etbs([(0, 2, 1.5, 2.), (1, 2, 1, 2.)])
[[0, [6.0, 1]], [0, [3.0, 1]], [1, [1., 1]], [1, [2., 1]]]
pyscf.gto.mole.expand_etb(l, n, alpha, beta)

Generate the exponents of even tempered basis for Mole.basis. .. math:

e = e^{-\alpha * \beta^{i-1}} for i = 1 .. n
Args:
lint

Angular momentum

nint

Number of GTOs

Returns:

Formated basis

Examples:

>>> gto.expand_etb(1, 3, 1.5, 2)
[[1, [6.0, 1]], [1, [3.0, 1]], [1, [1.5, 1]]]
pyscf.gto.mole.expand_etbs(etbs)

Generate even tempered basis. See also expand_etb()

Args:

etbs = [(l, n, alpha, beta), (l, n, alpha, beta),…]

Returns:

Formated basis

Examples:

>>> gto.expand_etbs([(0, 2, 1.5, 2.), (1, 2, 1, 2.)])
[[0, [6.0, 1]], [0, [3.0, 1]], [1, [1., 1]], [1, [2., 1]]]
pyscf.gto.mole.fakemol_for_charges(coords, expnt=1e+16)

Construct a fake Mole object that holds the charges on the given coordinates (coords). The shape of the charge can be a normal distribution with the Gaussian exponent (expnt).

pyscf.gto.mole.filatov_nuc_mod(nuc_charge, nucprop={})

Generate the nuclear charge distribution parameter zeta rho(r) = nuc_charge * Norm * exp(-zeta * r^2)

Ref. M. Filatov and D. Cremer, Theor. Chem. Acc. 108, 168 (2002)
  1. Filatov and D. Cremer, Chem. Phys. Lett. 351, 259 (2002)

pyscf.gto.mole.format_atom(atoms, origin=0, axes=None, unit='angstrom')

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 spearate 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 reserve 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]]]
pyscf.gto.mole.format_basis(basis_tab)

Convert the input Mole.basis to the internal data format.

``{ atom: [(l, ((-exp, c_1, c_2, ..),

(-exp, c_1, c_2, ..))),

(l, ((-exp, c_1, c_2, ..),

(-exp, c_1, c_2, ..)))], … }``

Args:
basis_tabdict

Similar to Mole.basis, it cannot be a str

Returns:

Formated basis

Examples:

>>> gto.format_basis({'H':'sto-3g', 'H^2': '3-21g'})
{'H': [[0,
    [3.4252509099999999, 0.15432897000000001],
    [0.62391373000000006, 0.53532813999999995],
    [0.16885539999999999, 0.44463454000000002]]],
 'H^2': [[0,
    [5.4471780000000001, 0.15628500000000001],
    [0.82454700000000003, 0.90469100000000002]],
    [0, [0.18319199999999999, 1.0]]]}
pyscf.gto.mole.format_ecp(ecp_tab)

Convert the input ecp (dict) to the internal data format:

{ atom: (nelec,  # core electrons
((l, # l=-1 for UL, l>=0 for Ul to indicate |l><l|
(((exp_1, c_1), # for r^0

(exp_2, c_2), …),

((exp_1, c_1), # for r^1

(exp_2, c_2), …),

((exp_1, c_1), # for r^2

…))))),

…}

pyscf.gto.mole.from_zmatrix(atomstr)
>>> a = """H
H 1 2.67247631453057
H 1 4.22555607338457 2 50.7684795164077
H 1 2.90305235726773 2 79.3904651036893 3 6.20854462618583"""
>>> for x in zmat2cart(a): print(x)
['H', array([ 0.,  0.,  0.])]
['H', array([ 2.67247631,  0.        ,  0.        ])]
['H', array([ 2.67247631,  0.        ,  3.27310166])]
['H', array([ 0.53449526,  0.30859098,  2.83668811])]
pyscf.gto.mole.fromfile(filename, format=None)

Read molecular geometry from a file (in testing)

Supported formats:
raw: Each line is <symobl> <x> <y> <z>
xyz: XYZ cartesian coordinates format
zmat: Z-matrix format
pyscf.gto.mole.fromstring(string, format='xyz')

Convert the string of the specified format to internal format (in testing)

Supported formats:
raw: Each line is <symobl> <x> <y> <z>
xyz: XYZ cartesian coordinates format
zmat: Z-matrix format
pyscf.gto.mole.gaussian_int(n, alpha)

int_0^inf x^n exp(-alpha x^2) dx

pyscf.gto.mole.gto_norm(l, expnt)

Normalized factor for GTO radial part \(g=r^l e^{-\alpha r^2}\)

\[\frac{1}{\sqrt{\int g^2 r^2 dr}} = \sqrt{\frac{2^{2l+3} (l+1)! (2a)^{l+1.5}}{(2l+2)!\sqrt{\pi}}}\]

Ref: H. B. Schlegel and M. J. Frisch, Int. J. Quant. Chem., 54(1995), 83-87.

Args:
l (int):

angular momentum

expnt :

exponent \(\alpha\)

Returns:

normalization factor

Examples:

>>> print(gto_norm(0, 1))
2.5264751109842591
pyscf.gto.mole.inertia_moment(mol, mass=None, coords=None)
pyscf.gto.mole.inter_distance(mol, coords=None)

Inter-particle distance array

pyscf.gto.mole.intor_cross(intor, mol1, mol2, comp=None)

1-electron integrals from two molecules like

\[\langle \mu | intor | \nu \rangle, \mu \in mol1, \nu \in mol2\]
Args:
intorstr

Name of the 1-electron integral, such as int1e_ovlp_sph (spherical overlap), int1e_nuc_cart (cartesian nuclear attraction), int1e_ipovlp_spinor (spinor overlap gradients), etc. Ref to getints() for the full list of available 1-electron integral names

mol1, mol2:

Mole objects

Kwargs:
compint

Components of the integrals, e.g. int1e_ipovlp_sph has 3 components

Returns:

ndarray of 1-electron integrals, can be either 2-dim or 3-dim, depending on comp

Examples:

Compute the overlap between H2 molecule and O atom

>>> mol1 = gto.M(atom='H 0 1 0; H 0 0 1', basis='sto3g')
>>> mol2 = gto.M(atom='O 0 0 0', basis='sto3g')
>>> gto.intor_cross('int1e_ovlp_sph', mol1, mol2)
[[ 0.04875181  0.44714688  0.          0.37820346  0.        ]
 [ 0.04875181  0.44714688  0.          0.          0.37820346]]
pyscf.gto.mole.is_same_mol(mol1, mol2, tol=1e-05, cmp_basis=True, ignore_chiral=False)

Compare the two molecules whether they have the same structure.

Kwargs:
tolfloat

In Bohr

cmp_basisbool

Whether to compare basis functions for the two molecules

pyscf.gto.mole.len_cart(l)

The number of Cartesian function associated with given angular momentum.

pyscf.gto.mole.len_spinor(l, kappa)

The number of spinor associated with given angular momentum and kappa. If kappa is 0, return 4l+2

pyscf.gto.mole.loads(molstr)

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

pyscf.gto.mole.make_atm_env(atom, ptr=0, nuclear_model=1, nucprop={})

Convert the internal format Mole._atom to the format required by libcint integrals

pyscf.gto.mole.make_bas_env(basis_add, atom_id=0, ptr=0)

Convert Mole.basis to the argument bas for libcint integrals

pyscf.gto.mole.make_ecp_env(mol, _atm, ecp, pre_env=[])

Generate the input arguments _ecpbas for ECP integrals

pyscf.gto.mole.make_env(atoms, basis, pre_env=[], nucmod={}, nucprop={})

Generate the input arguments for libcint library based on internal format Mole._atom and Mole._basis

pyscf.gto.mole.nao_2c(mol)

Total number of contracted spinor GTOs for the given Mole object

pyscf.gto.mole.nao_2c_range(mol, bas_id0, bas_id1)

Lower and upper boundary of contracted spinor basis functions associated with the given shell range

Args:
mol :

Mole object

bas_id0int

start shell id, 0-based

bas_id1int

stop shell id, 0-based

Returns:

tupel of start basis function id and the stop function id

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.nao_2c_range(mol, 2, 4)
(4, 12)
pyscf.gto.mole.nao_cart(mol)

Total number of contracted cartesian GTOs for the given Mole object

pyscf.gto.mole.nao_nr(mol, cart=None)

Total number of contracted GTOs for the given Mole object

pyscf.gto.mole.nao_nr_range(mol, bas_id0, bas_id1)

Lower and upper boundary of contracted spherical basis functions associated with the given shell range

Args:
mol :

Mole object

bas_id0int

start shell id

bas_id1int

stop shell id

Returns:

tupel of start basis function id and the stop function id

Examples:

>>> mol = gto.M(atom='O 0 0 0; C 0 0 1', basis='6-31g')
>>> gto.nao_nr_range(mol, 2, 4)
(2, 6)
pyscf.gto.mole.npgto_nr(mol, cart=None)

Total number of primitive spherical GTOs for the given Mole object

pyscf.gto.mole.offset_2c_by_atom(mol)

2-component AO offset for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

pyscf.gto.mole.offset_nr_by_atom(mol, ao_loc=None)

AO offsets for each atom. Return a list, each item of the list gives (start-shell-id, stop-shell-id, start-AO-id, stop-AO-id)

pyscf.gto.mole.pack(mol)

Pack the input args of Mole to a dict.

Note this function only pack the input arguments (not the entire Mole class). Modifications to mol._atm, mol._bas, mol._env are not tracked. Use dumps() to serialize the entire Mole object.

pyscf.gto.mole.same_basis_set(mol1, mol2)

Check whether two molecules use the same basis sets. The two molecules can have different geometry.

pyscf.gto.mole.same_mol(mol1, mol2, tol=1e-05, cmp_basis=True, ignore_chiral=False)

Compare the two molecules whether they have the same structure.

Kwargs:
tolfloat

In Bohr

cmp_basisbool

Whether to compare basis functions for the two molecules

pyscf.gto.mole.search_ao_label(mol, label)

Find the index of the AO basis function based on the given ao_label

Args:
ao_labelstring or a list of strings

The regular expression pattern to match the orbital labels returned by mol.ao_labels()

Returns:

A list of index for the AOs that matches the given ao_label RE pattern

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='ccpvtz')
>>> mol.search_ao_label('Cl.*p')
[19 20 21 22 23 24 25 26 27 28 29 30]
>>> mol.search_ao_label('Cl 2p')
[19 20 21]
>>> mol.search_ao_label(['Cl.*d', 'Cl 4p'])
[25 26 27 31 32 33 34 35 36 37 38 39 40]
pyscf.gto.mole.search_ao_nr(mol, atm_id, l, m, atmshell)

Search the first basis function id (not the shell id) which matches the given atom-id, angular momentum magnetic angular momentum, principal shell.

Args:
atm_idint

atom id, 0-based

lint

angular momentum

mint

magnetic angular momentum

atmshellint

principal quantum number

Returns:

basis function id, 0-based. If not found, return None

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='sto-3g')
>>> mol.search_ao_nr(1, 1, -1, 3) # Cl 3px
7
pyscf.gto.mole.search_ao_r(mol, atm_id, l, j, m, atmshell)
pyscf.gto.mole.search_shell_id(mol, atm_id, l)

Search the first basis/shell id (not the basis function id) which matches the given atom-id and angular momentum

Args:
atm_idint

atom id, 0-based

lint

angular momentum

Returns:

basis id, 0-based. If not found, return None

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='sto-3g')
>>> mol.search_shell_id(1, 1) # Cl p shell
4
>>> mol.search_shell_id(1, 2) # Cl d shell
None
pyscf.gto.mole.sph2spinor_kappa(kappa, l=None)

Real spherical to spinor transformation matrix for kappa

pyscf.gto.mole.sph2spinor_l(l)

Real spherical to spinor transformation matrix for angular moment l

pyscf.gto.mole.sph_labels(mol, fmt=True, base=0)

Labels for spherical GTO functions

Kwargs:

fmt : str or bool if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-real-spherical-notation] or formatted strings based on the argument “fmt”

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='sto-3g')
>>> gto.sph_labels(mol)
[(0, 'H', '1s', ''), (1, 'Cl', '1s', ''), (1, 'Cl', '2s', ''), (1, 'Cl', '3s', ''),
 (1, 'Cl', '2p', 'x'), (1, 'Cl', '2p', 'y'), (1, 'Cl', '2p', 'z'), (1, 'Cl', '3p', 'x'),
 (1, 'Cl', '3p', 'y'), (1, 'Cl', '3p', 'z')]
pyscf.gto.mole.spheric_labels(mol, fmt=True, base=0)

Labels for spherical GTO functions

Kwargs:

fmt : str or bool if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-real-spherical-notation] or formatted strings based on the argument “fmt”

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='sto-3g')
>>> gto.sph_labels(mol)
[(0, 'H', '1s', ''), (1, 'Cl', '1s', ''), (1, 'Cl', '2s', ''), (1, 'Cl', '3s', ''),
 (1, 'Cl', '2p', 'x'), (1, 'Cl', '2p', 'y'), (1, 'Cl', '2p', 'z'), (1, 'Cl', '3p', 'x'),
 (1, 'Cl', '3p', 'y'), (1, 'Cl', '3p', 'z')]
pyscf.gto.mole.spherical_labels(mol, fmt=True, base=0)

Labels for spherical GTO functions

Kwargs:

fmt : str or bool if fmt is boolean, it controls whether to format the labels and the default format is “%d%3s %s%-4s”. if fmt is string, the string will be used as the print format.

Returns:

List of [(atom-id, symbol-str, nl-str, str-of-real-spherical-notation] or formatted strings based on the argument “fmt”

Examples:

>>> mol = gto.M(atom='H 0 0 0; Cl 0 0 1', basis='sto-3g')
>>> gto.sph_labels(mol)
[(0, 'H', '1s', ''), (1, 'Cl', '1s', ''), (1, 'Cl', '2s', ''), (1, 'Cl', '3s', ''),
 (1, 'Cl', '2p', 'x'), (1, 'Cl', '2p', 'y'), (1, 'Cl', '2p', 'z'), (1, 'Cl', '3p', 'x'),
 (1, 'Cl', '3p', 'y'), (1, 'Cl', '3p', 'z')]
pyscf.gto.mole.spinor_labels(mol, fmt=True, base=0)

Labels of spinor GTO functions

pyscf.gto.mole.time_reversal_map(mol)

The index to map the spinor functions and its time reversal counterpart. The returned indices have postive or negative values. For the i-th basis function, if the returned j = idx[i] < 0, it means \(T|i\rangle = -|j\rangle\), otherwise \(T|i\rangle = |j\rangle\)

pyscf.gto.mole.to_uncontracted_cartesian_basis(mol)

Decontract the basis of a Mole or a Cell. Returns a Mole (Cell) object with the uncontracted basis environment and a list of coefficients that transform the uncontracted cartesian basis to the original basis. Each element in the list corresponds to one shell of the original Mole (Cell).

Examples:

>>> mol = gto.M(atom='Ne', basis='ccpvdz')
>>> pmol, ctr_coeff = mol.to_uncontracted_cartesian_basis()
>>> c = scipy.linalg.block_diag(*ctr_coeff)
>>> s = reduce(numpy.dot, (c.T, pmol.intor('int1e_ovlp'), c))
>>> abs(s-mol.intor('int1e_ovlp')).max()
0.0
pyscf.gto.mole.tofile(mol, filename, format=None)

Write molecular geometry to a file of the required format.

Supported output formats:
raw: Each line is <symobl> <x> <y> <z>
xyz: XYZ cartesian coordinates format
zmat: Z-matrix format
pyscf.gto.mole.tostring(mol, format='raw')

Convert molecular geometry to a string of the required format.

Supported output formats:
raw: Each line is <symobl> <x> <y> <z>
xyz: XYZ cartesian coordinates format
zmat: Z-matrix format
pyscf.gto.mole.tot_electrons(mol)

Total number of electrons for the given molecule

Returns:

electron number in integer

Examples:

>>> mol = gto.M(atom='H 0 1 0; C 0 0 1', charge=1)
>>> mol.tot_electrons()
6
pyscf.gto.mole.uncontract(_basis)

Uncontract internal format _basis

Examples:

>>> gto.uncontract(gto.load('sto3g', 'He'))
[[0, [6.36242139, 1]], [0, [1.158923, 1]], [0, [0.31364979, 1]]]
pyscf.gto.mole.uncontracted_basis(_basis)

Uncontract internal format _basis

Examples:

>>> gto.uncontract(gto.load('sto3g', 'He'))
[[0, [6.36242139, 1]], [0, [1.158923, 1]], [0, [0.31364979, 1]]]
pyscf.gto.mole.unpack(moldic)

Unpack a dict which is packed by pack(), to generate the input arguments for Mole object.

pyscf.gto.mole.zmat(atomstr)
>>> a = """H
H 1 2.67247631453057
H 1 4.22555607338457 2 50.7684795164077
H 1 2.90305235726773 2 79.3904651036893 3 6.20854462618583"""
>>> for x in zmat2cart(a): print(x)
['H', array([ 0.,  0.,  0.])]
['H', array([ 2.67247631,  0.        ,  0.        ])]
['H', array([ 2.67247631,  0.        ,  3.27310166])]
['H', array([ 0.53449526,  0.30859098,  2.83668811])]
pyscf.gto.mole.zmat2cart(atomstr)
>>> a = """H
H 1 2.67247631453057
H 1 4.22555607338457 2 50.7684795164077
H 1 2.90305235726773 2 79.3904651036893 3 6.20854462618583"""
>>> for x in zmat2cart(a): print(x)
['H', array([ 0.,  0.,  0.])]
['H', array([ 2.67247631,  0.        ,  0.        ])]
['H', array([ 2.67247631,  0.        ,  3.27310166])]
['H', array([ 0.53449526,  0.30859098,  2.83668811])]

pyscf.gto.moleintor module

A low level interface to libcint library. It’s recommended to use the Mole.intor method to drive the integral evaluation funcitons.

pyscf.gto.moleintor.ascint3(intor_name)

convert cint2 function name to cint3 function name

pyscf.gto.moleintor.getints(intor_name, atm, bas, env, shls_slice=None, comp=None, hermi=0, aosym='s1', ao_loc=None, cintopt=None, out=None)

1e and 2e integral generator.

Args:

intor_name : str

Function

Expression

“int1e_ovlp”

( | )

“int1e_nuc”

( | nuc | )

“int1e_kin”

(.5 | p dot p)

“int1e_ia01p”

(#C(0 1) | nabla-rinv | cross p)

“int1e_giao_irjxp”

(#C(0 1) | r cross p)

“int1e_cg_irxp”

(#C(0 1) | rc cross p)

“int1e_giao_a11part”

(-.5 | nabla-rinv | r)

“int1e_cg_a11part”

(-.5 | nabla-rinv | rc)

“int1e_a01gp”

(g | nabla-rinv cross p |)

“int1e_igkin”

(#C(0 .5) g | p dot p)

“int1e_igovlp”

(#C(0 1) g |)

“int1e_ignuc”

(#C(0 1) g | nuc |)

“int1e_z”

( | zc | )

“int1e_zz”

( | zc zc | )

“int1e_r”

( | rc | )

“int1e_r2”

( | rc dot rc | )

“int1e_rr”

( | rc rc | )

“int1e_rrr”

( | rc rc rc | )

“int1e_rrrr”

( | rc rc rc rc | )

“int1e_pnucp”

(p* | nuc dot p | )

“int1e_prinvxp”

(p* | rinv cross p | )

“int1e_ipovlp”

(nabla |)

“int1e_ipkin”

(.5 nabla | p dot p)

“int1e_ipnuc”

(nabla | nuc |)

“int1e_iprinv”

(nabla | rinv |)

“int1e_rinv”

(| rinv |)

“int1e_pnucxp”

(p* | nuc cross p | )

“int1e_irp”

( | rc nabla | )

“int1e_irrp”

( | rc rc nabla | )

“int1e_irpr”

( | rc nabla rc | )

“int1e_ggovlp”

( | g g | )

“int1e_ggkin”

(.5 | g g p dot p | )

“int1e_ggnuc”

( | g g nuc | )

“int1e_grjxp”

( | g r cross p | )

“ECPscalar”

AREP ECP integrals, similar to int1e_nuc

“ECPscalar_ipnuc”

(nabla i | ECP | ), similar to int1e_ipnuc

“ECPscalar_iprinv”

similar to int1e_iprinv for a specific atom

“ECPscalar_ignuc”

similar to int1e_ignuc

“ECPscalar_iprinvip”

similar to int1e_iprinvip

“ECPso”

< | Spin-orbit ECP | >

“int1e_ovlp_spinor”

( | )

“int1e_nuc_spinor”

( | nuc |)

“int1e_srsr_spinor”

(sigma dot r | sigma dot r)

“int1e_sr_spinor”

(sigma dot r |)

“int1e_srsp_spinor”

(sigma dot r | sigma dot p)

“int1e_spsp_spinor”

(sigma dot p | sigma dot p)

“int1e_sp_spinor”

(sigma dot p |)

“int1e_spnucsp_spinor”

(sigma dot p | nuc | sigma dot p)

“int1e_srnucsr_spinor”

(sigma dot r | nuc | sigma dot r)

“int1e_govlp_spinor”

(g |)

“int1e_gnuc_spinor”

(g | nuc |)

“int1e_cg_sa10sa01_spinor”

(.5 sigma cross rc | sigma cross nabla-rinv |)

“int1e_cg_sa10sp_spinor”

(.5 rc cross sigma | sigma dot p)

“int1e_cg_sa10nucsp_spinor”

(.5 rc cross sigma | nuc | sigma dot p)

“int1e_giao_sa10sa01_spinor”

(.5 sigma cross r | sigma cross nabla-rinv |)

“int1e_giao_sa10sp_spinor”

(.5 r cross sigma | sigma dot p)

“int1e_giao_sa10nucsp_spinor”

(.5 r cross sigma | nuc | sigma dot p)

“int1e_sa01sp_spinor”

(| nabla-rinv cross sigma | sigma dot p)

“int1e_spgsp_spinor”

(g sigma dot p | sigma dot p)

“int1e_spgnucsp_spinor”

(g sigma dot p | nuc | sigma dot p)

“int1e_spgsa01_spinor”

(g sigma dot p | nabla-rinv cross sigma |)

“int1e_spspsp_spinor”

(sigma dot p | sigma dot p sigma dot p)

“int1e_spnuc_spinor”

(sigma dot p | nuc |)

“int1e_ipovlp_spinor”

(nabla |)

“int1e_ipkin_spinor”

(.5 nabla | p dot p)

“int1e_ipnuc_spinor”

(nabla | nuc |)

“int1e_iprinv_spinor”

(nabla | rinv |)

“int1e_ipspnucsp_spinor”

(nabla sigma dot p | nuc | sigma dot p)

“int1e_ipsprinvsp_spinor”

(nabla sigma dot p | rinv | sigma dot p)

“int2e”

( , | , )

“int2e_ig1”

(#C(0 1) g , | , )

“int2e_gg1”

(g g , | , )

“int2e_g1g2”

(g , | g , )

“int2e_p1vxp1”

( p* , cross p | , ) ; SSO

“int2e_spinor”

(, | , )

“int2e_spsp1_spinor”

(sigma dot p , sigma dot p | , )

“int2e_spsp1spsp2_spinor”

(sigma dot p , sigma dot p | sigma dot p , sigma dot p )

“int2e_srsr1_spinor”

(sigma dot r , sigma dot r | ,)

“int2e_srsr1srsr2_spinor”

(sigma dot r , sigma dot r | sigma dot r , sigma dot r)

“int2e_cg_sa10sp1_spinor”

(.5 rc cross sigma , sigma dot p | ,)

“int2e_cg_sa10sp1spsp2_spinor”

(.5 rc cross sigma , sigma dot p | sigma dot p , sigma dot p )

“int2e_giao_sa10sp1_spinor”

(.5 r cross sigma , sigma dot p | ,)

“int2e_giao_sa10sp1spsp2_spinor”

(.5 r cross sigma , sigma dot p | sigma dot p , sigma dot p )

“int2e_g1_spinor”

(g , | ,)

“int2e_spgsp1_spinor”

(g sigma dot p , sigma dot p | ,)

“int2e_g1spsp2_spinor”

(g , | sigma dot p , sigma dot p)

“int2e_spgsp1spsp2_spinor”

(g sigma dot p , sigma dot p | sigma dot p , sigma dot p)

“int2e_spv1_spinor”

(sigma dot p , | ,)

“int2e_vsp1_spinor”

(, sigma dot p | ,)

“int2e_spsp2_spinor”

(, | sigma dot p , sigma dot p)

“int2e_spv1spv2_spinor”

(sigma dot p , | sigma dot p ,)

“int2e_vsp1spv2_spinor”

(, sigma dot p | sigma dot p ,)

“int2e_spv1vsp2_spinor”

(sigma dot p , | , sigma dot p)

“int2e_vsp1vsp2_spinor”

(, sigma dot p | , sigma dot p)

“int2e_spv1spsp2_spinor”

(sigma dot p , | sigma dot p , sigma dot p)

“int2e_vsp1spsp2_spinor”

(, sigma dot p | sigma dot p , sigma dot p)

“int2e_ig1”

(#C(0 1) g , | , )

“int2e_ip1”

(nabla , | ,)

“int2e_ip1_spinor”

(nabla , | ,)

“int2e_ipspsp1_spinor”

(nabla sigma dot p , sigma dot p | ,)

“int2e_ip1spsp2_spinor”

(nabla , | sigma dot p , sigma dot p)

“int2e_ipspsp1spsp2_spinor”

(nabla sigma dot p , sigma dot p | sigma dot p , sigma dot p)

“int2e_ipsrsr1_spinor”

(nabla sigma dot r , sigma dot r | ,)

“int2e_ip1srsr2_spinor”

(nabla , | sigma dot r , sigma dot r)

“int2e_ipsrsr1srsr2_spinor”

(nabla sigma dot r , sigma dot r | sigma dot r , sigma dot r)

“int2e_ip1”

(nabla , | ,)

“int2e_ssp1ssp2_spinor”

( , sigma dot p | gaunt | , sigma dot p)

“int2e_cg_ssa10ssp2_spinor”

(rc cross sigma , | gaunt | , sigma dot p)

“int2e_giao_ssa10ssp2_spinor”

(r cross sigma , | gaunt | , sigma dot p)

“int2e_gssp1ssp2_spinor”

(g , sigma dot p | gaunt | , sigma dot p)

“int2e_ipip1”

( nabla nabla , | , )

“int2e_ipvip1”

( nabla , nabla | , )

“int2e_ip1ip2”

( nabla , | nabla , )

“int3c2e_ip1”

(nabla , | )

“int3c2e_ip2”

( , | nabla)

“int2c2e_ip1”

(nabla | r12 | )

“int3c2e_spinor”

(nabla , | )

“int3c2e_spsp1_spinor”

(nabla , | )

“int3c2e_ip1_spinor”

(nabla , | )

“int3c2e_ip2_spinor”

( , | nabla)

“int3c2e_ipspsp1_spinor”

(nabla sigma dot p , sigma dot p | )

“int3c2e_spsp1ip2_spinor”

(sigma dot p , sigma dot p | nabla )

“ECPscalar_spinor”

AREP ECP integrals, similar to int1e_nuc

“ECPscalar_ipnuc_spinor”

(nabla i | ECP | ), similar to int1e_ipnuc

“ECPscalar_iprinv_spinor”

similar to int1e_iprinv for a specific atom

“ECPscalar_ignuc_spinor”

similar to int1e_ignuc

“ECPscalar_iprinvip_spinor”

similar to int1e_iprinvip

“ECPso_spinor”

< | sigam dot Spin-orbit ECP | >

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

Kwargs:
shls_slice8-element list

(ish_start, ish_end, jsh_start, jsh_end, ksh_start, ksh_end, lsh_start, lsh_end)

compint

Components of the integrals, e.g. int1e_ipovlp has 3 components.

hermiint (1e integral only)

Symmetry of the 1e integrals

0 : no symmetry assumed (default)
1 : hermitian
2 : anti-hermitian
aosymstr (2e integral only)

Symmetry of the 2e integrals

4 or ‘4’ or ‘s4’: 4-fold symmetry (default)
‘2ij’ or ‘s2ij’ : symmetry between i, j in (ij|kl)
‘2kl’ or ‘s2kl’ : symmetry between k, l in (ij|kl)
1 or ‘1’ or ‘s1’: no symmetry
outndarray (2e integral only)

array to store the 2e AO integrals

Returns:

ndarray of 1-electron integrals, can be either 2-dim or 3-dim, depending on comp

Examples:

>>> mol.build(atom='H 0 0 0; H 0 0 1.1', basis='sto-3g')
>>> gto.getints('int1e_ipnuc_sph', mol._atm, mol._bas, mol._env, comp=3) # <nabla i | V_nuc | j>
[[[ 0.          0.        ]
  [ 0.          0.        ]]
 [[ 0.          0.        ]
  [ 0.          0.        ]]
 [[ 0.10289944  0.48176097]
  [-0.48176097 -0.10289944]]]
pyscf.gto.moleintor.getints2c(intor_name, atm, bas, env, shls_slice=None, comp=1, hermi=0, ao_loc=None, cintopt=None, out=None)
pyscf.gto.moleintor.getints3c(intor_name, atm, bas, env, shls_slice=None, comp=1, aosym='s1', ao_loc=None, cintopt=None, out=None)
pyscf.gto.moleintor.getints4c(intor_name, atm, bas, env, shls_slice=None, comp=1, aosym='s1', ao_loc=None, cintopt=None, out=None)
pyscf.gto.moleintor.getints_by_shell(intor_name, shls, atm, bas, env, comp=1)

For given 2, 3 or 4 shells, interface for libcint to get 1e, 2e, 2-center-2e or 3-center-2e integrals

Args:
intor_namestr

See also getints() for the supported intor_name

shlslist of int

The AO shell-ids of the integrals

atmint32 ndarray

libcint integral function argument

basint32 ndarray

libcint integral function argument

envfloat64 ndarray

libcint integral function argument

Kwargs:
compint

Components of the integrals, e.g. int1e_ipovlp has 3 components.

Returns:

ndarray of 2-dim to 5-dim, depending on the integral type (1e, 2e, 3c-2e, 2c2e) and the value of comp

Examples:

The gradients of the spherical 2e integrals

>>> mol.build(atom='H 0 0 0; H 0 0 1.1', basis='sto-3g')
>>> gto.getints_by_shell('int2e_ip1_sph', (0,1,0,1), mol._atm, mol._bas, mol._env, comp=3)
[[[[[-0.        ]]]]
  [[[[-0.        ]]]]
  [[[[-0.08760462]]]]]
pyscf.gto.moleintor.make_cintopt(atm, bas, env, intor)
pyscf.gto.moleintor.make_loc(bas, key)

pyscf.gto.pp_int module

Analytic GTH-PP integrals for open boundary conditions.

See also pyscf/pbc/gto/pseudo/pp_int.py

pyscf.gto.pp_int.get_gth_pp(mol)

Module contents