pyscf.ao2mo package¶
Submodules¶
pyscf.ao2mo.addons module¶
-
class
pyscf.ao2mo.addons.
load
(eri, dataname='eri_mo')¶ Bases:
object
load 2e integrals from hdf5 file
- Usage:
- with load(erifile) as eri:
print(eri.shape)
-
pyscf.ao2mo.addons.
restore
(symmetry, eri, norb, tao=None)¶ Convert the 2e integrals (in Chemist’s notation) between different level of permutation symmetry (8-fold, 4-fold, or no symmetry)
- Args:
- symmetryint or str
code to present the target symmetry of 2e integrals
‘s8’ or ‘8’ or 8 : 8-fold symmetry‘s4’ or ‘4’ or 4 : 4-fold symmetry‘s1’ or ‘1’ or 1 : no symmetry‘s2ij’ or ‘2ij’ : symmetric ij pair for (ij|kl) (TODO)‘s2ij’ or ‘2kl’ : symmetric kl pair for (ij|kl) (TODO)Note the 4-fold symmetry requires (ij|kl) == (ij|lk) == (ij|lk) while (ij|kl) != (kl|ij) is not required.
- erindarray
The symmetry of eri is determined by the size of eri and norb
- norbint
The symmetry of eri is determined by the size of eri and norb
- Returns:
ndarray. The shape depends on the target symmetry.
8 : (norb*(norb+1)/2)*(norb*(norb+1)/2+1)/24 : (norb*(norb+1)/2, norb*(norb+1)/2)1 : (norb, norb, norb, norb)
Examples:
>>> from pyscf import gto >>> from pyscf.scf import _vhf >>> from pyscf import ao2mo >>> mol = gto.M(atom='O 0 0 0; H 0 1 0; H 0 0 1', basis='sto3g') >>> eri = mol.intor('int2e') >>> eri1 = ao2mo.restore(1, eri, mol.nao_nr()) >>> eri4 = ao2mo.restore(4, eri, mol.nao_nr()) >>> eri8 = ao2mo.restore(8, eri, mol.nao_nr()) >>> print(eri1.shape) (7, 7, 7, 7) >>> print(eri1.shape) (28, 28) >>> print(eri1.shape) (406,)
pyscf.ao2mo.incore module¶
-
pyscf.ao2mo.incore.
full
(eri_ao, mo_coeff, verbose=0, compact=True, **kwargs)¶ MO integral transformation for the given orbital.
- Args:
- eri_aondarray
AO integrals, can be either 8-fold or 4-fold symmetry.
- mo_coeffndarray
Transform (ij|kl) with the same set of orbitals.
- Kwargs:
- verboseint
Print level
- compactbool
When compact is True, the returned MO integrals have 4-fold symmetry. Otherwise, return the “plain” MO integrals.
- Returns:
2D array of transformed MO integrals. The MO integrals may or may not have the permutation symmetry (controlled by the kwargs compact)
Examples:
>>> from pyscf import gto >>> from pyscf.scf import _vhf >>> from pyscf import ao2mo >>> mol = gto.M(atom='O 0 0 0; H 0 1 0; H 0 0 1', basis='sto3g') >>> eri = mol.intor('int2e_sph', aosym='s8') >>> mo1 = numpy.random.random((mol.nao_nr(), 10)) >>> eri1 = ao2mo.incore.full(eri, mo1) >>> print(eri1.shape) (55, 55) >>> eri1 = ao2mo.incore.full(eri, mo1, compact=False) >>> print(eri1.shape) (100, 100)
-
pyscf.ao2mo.incore.
general
(eri_ao, mo_coeffs, verbose=0, compact=True, **kwargs)¶ For the given four sets of orbitals, transfer the 8-fold or 4-fold 2e AO integrals to MO integrals.
- Args:
- eri_aondarray
AO integrals, can be either 8-fold or 4-fold symmetry.
- mo_coeffs4-item list of ndarray
Four sets of orbital coefficients, corresponding to the four indices of (ij|kl)
- Kwargs:
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
2D array of transformed MO integrals. The MO integrals may or may not have the permutation symmetry, depending on the given orbitals, and the kwargs compact. If the four sets of orbitals are identical, the MO integrals will at most have 4-fold symmetry.
Examples:
>>> from pyscf import gto >>> from pyscf.scf import _vhf >>> from pyscf import ao2mo >>> mol = gto.M(atom='O 0 0 0; H 0 1 0; H 0 0 1', basis='sto3g') >>> eri = mol.intor('int2e_sph', aosym='s8') >>> mo1 = numpy.random.random((mol.nao_nr(), 10)) >>> mo2 = numpy.random.random((mol.nao_nr(), 8)) >>> mo3 = numpy.random.random((mol.nao_nr(), 6)) >>> mo4 = numpy.random.random((mol.nao_nr(), 4)) >>> eri1 = ao2mo.incore.general(eri, (mo1,mo2,mo3,mo4)) >>> print(eri1.shape) (80, 24) >>> eri1 = ao2mo.incore.general(eri, (mo1,mo2,mo3,mo3)) >>> print(eri1.shape) (80, 21) >>> eri1 = ao2mo.incore.general(eri, (mo1,mo2,mo3,mo3), compact=False) >>> print(eri1.shape) (80, 36) >>> eri1 = ao2mo.incore.general(eri, (mo1,mo1,mo2,mo2)) >>> print(eri1.shape) (55, 36) >>> eri1 = ao2mo.incore.general(eri, (mo1,mo2,mo1,mo2)) >>> print(eri1.shape) (80, 80)
-
pyscf.ao2mo.incore.
half_e1
(eri_ao, mo_coeffs, compact=True)¶ Given two set of orbitals, half transform the (ij| pair of 8-fold or 4-fold AO integrals (ij|kl)
- Args:
- eri_aondarray
AO integrals, can be either 8-fold or 4-fold symmetry.
- mo_coeffslist of ndarray
Two sets of orbital coefficients, corresponding to the i, j indices of (ij|kl)
- Kwargs:
- compactbool
When compact is True, the returned MO integrals uses the highest possible permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
ndarray of transformed MO integrals. The MO integrals may or may not have the permutation symmetry, depending on the given orbitals, and the kwargs compact.
Examples:
>>> from pyscf import gto >>> from pyscf import ao2mo >>> mol = gto.M(atom='O 0 0 0; H 0 1 0; H 0 0 1', basis='sto3g') >>> eri = mol.intor('int2e_sph', aosym='s8') >>> mo1 = numpy.random.random((mol.nao_nr(), 10)) >>> mo2 = numpy.random.random((mol.nao_nr(), 8)) >>> eri1 = ao2mo.incore.half_e1(eri, (mo1,mo2)) >>> eri1 = ao2mo.incore.half_e1(eri, (mo1,mo2)) >>> print(eri1.shape) (80, 28) >>> eri1 = ao2mo.incore.half_e1(eri, (mo1,mo2), compact=False) >>> print(eri1.shape) (80, 28) >>> eri1 = ao2mo.incore.half_e1(eri, (mo1,mo1)) >>> print(eri1.shape) (55, 28)
-
pyscf.ao2mo.incore.
iden_coeffs
(mo1, mo2)¶
pyscf.ao2mo.outcore module¶
-
pyscf.ao2mo.outcore.
balance_partition
(ao_loc, blksize, start_id=0, stop_id=None)¶
-
pyscf.ao2mo.outcore.
balance_segs
(segs_lst, blksize, start_id=0, stop_id=None)¶
-
pyscf.ao2mo.outcore.
full
(mol, mo_coeff, erifile, dataname='eri_mo', intor='int2e', aosym='s4', comp=None, max_memory=2000, ioblk_size=256, verbose=2, compact=True)¶ Transfer arbitrary spherical AO integrals to MO integrals for given orbitals
- Args:
- mol
Mole
object AO integrals will be generated in terms of mol._atm, mol._bas, mol._env
- mo_coeffndarray
Transform (ij|kl) with the same set of orbitals.
- erifilestr or h5py File or h5py Group object
To store the transformed integrals, in HDF5 format.
- mol
- Kwargs:
- datanamestr
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 dataname, the new integrals data will overwrite the old one.
- intorstr
Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Components of the integrals, e.g. int2e_ip_sph has 3 components.
- max_memoryfloat or int
The maximum size of cache to use (in MB), large cache may not improve performance.
- ioblk_sizefloat or int
The block size for IO, large block size may not improve performance
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
None
Examples:
>>> from pyscf import gto >>> from pyscf import ao2mo >>> import h5py >>> def view(h5file, dataname='eri_mo'): ... f5 = h5py.File(h5file, 'r') ... print('dataset %s, shape %s' % (str(f5.keys()), str(f5[dataname].shape))) ... f5.close() >>> mol = gto.M(atom='O 0 0 0; H 0 1 0; H 0 0 1', basis='sto3g') >>> mo1 = numpy.random.random((mol.nao_nr(), 10)) >>> ao2mo.outcore.full(mol, mo1, 'full.h5') >>> view('full.h5') dataset ['eri_mo'], shape (55, 55) >>> ao2mo.outcore.full(mol, mo1, 'full.h5', dataname='new', compact=False) >>> view('full.h5', 'new') dataset ['eri_mo', 'new'], shape (100, 100) >>> ao2mo.outcore.full(mol, mo1, 'full.h5', intor='int2e_ip1_sph', aosym='s1', comp=3) >>> view('full.h5') dataset ['eri_mo', 'new'], shape (3, 100, 100) >>> ao2mo.outcore.full(mol, mo1, 'full.h5', intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> view('full.h5') dataset ['eri_mo', 'new'], shape (3, 100, 55)
-
pyscf.ao2mo.outcore.
full_iofree
(mol, mo_coeff, intor='int2e', aosym='s4', comp=None, max_memory=2000, ioblk_size=256, verbose=2, compact=True)¶ Transfer arbitrary spherical AO integrals to MO integrals for given orbitals This function is a wrap for
ao2mo.outcore.general()
. It’s not really IO free. The returned MO integrals are held in memory. For backward compatibility, it is used to replace the non-existed function direct.full_iofree.- Args:
- mol
Mole
object AO integrals will be generated in terms of mol._atm, mol._bas, mol._env
- mo_coeffndarray
Transform (ij|kl) with the same set of orbitals.
- erifilestr
To store the transformed integrals, in HDF5 format.
- mol
- Kwargs
- datanamestr
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 dataname, the new integrals data will overwrite the old one.
- intorstr
Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Components of the integrals, e.g. int2e_ip_sph has 3 components.
- verboseint
Print level
- max_memoryfloat or int
The maximum size of cache to use (in MB), large cache may not improve performance.
- ioblk_sizefloat or int
The block size for IO, large block size may not improve performance
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
2D/3D MO-integral array. They may or may not have the permutation symmetry, depending on the given orbitals, and the kwargs compact. If the four sets of orbitals are identical, the MO integrals will at most have 4-fold symmetry.
Examples:
>>> from pyscf import gto >>> from pyscf import ao2mo >>> mol = gto.M(atom='O 0 0 0; H 0 1 0; H 0 0 1', basis='sto3g') >>> mo1 = numpy.random.random((mol.nao_nr(), 10)) >>> eri1 = ao2mo.outcore.full_iofree(mol, mo1) >>> print(eri1.shape) (55, 55) >>> eri1 = ao2mo.outcore.full_iofree(mol, mo1, compact=False) >>> print(eri1.shape) (100, 100) >>> eri1 = ao2mo.outcore.full_iofree(mol, mo1, intor='int2e_ip1_sph', aosym='s1', comp=3) >>> print(eri1.shape) (3, 100, 100) >>> eri1 = ao2mo.outcore.full_iofree(mol, mo1, intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> print(eri1.shape) (3, 100, 55)
-
pyscf.ao2mo.outcore.
general
(mol, mo_coeffs, erifile, dataname='eri_mo', intor='int2e', aosym='s4', comp=None, max_memory=2000, ioblk_size=256, verbose=2, compact=True)¶ For the given four sets of orbitals, transfer arbitrary spherical AO integrals to MO integrals on the fly.
- Args:
- mol
Mole
object AO integrals will be generated in terms of mol._atm, mol._bas, mol._env
- mo_coeffs4-item list of ndarray
Four sets of orbital coefficients, corresponding to the four indices of (ij|kl)
- erifilestr or h5py File or h5py Group object
To store the transformed integrals, in HDF5 format.
- mol
- Kwargs
- datanamestr
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 dataname, the new integrals data will overwrite the old one.
- intorstr
Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Components of the integrals, e.g. int2e_ip_sph has 3 components.
- max_memoryfloat or int
The maximum size of cache to use (in MB), large cache may not improve performance.
- ioblk_sizefloat or int
The block size for IO, large block size may not improve performance
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
None
Examples:
>>> from pyscf import gto >>> from pyscf import ao2mo >>> import h5py >>> def view(h5file, dataname='eri_mo'): ... f5 = h5py.File(h5file, 'r') ... print('dataset %s, shape %s' % (str(f5.keys()), str(f5[dataname].shape))) ... f5.close() >>> mol = gto.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)) >>> mo3 = numpy.random.random((mol.nao_nr(), 6)) >>> mo4 = numpy.random.random((mol.nao_nr(), 4)) >>> ao2mo.outcore.general(mol, (mo1,mo2,mo3,mo4), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 24) >>> ao2mo.outcore.general(mol, (mo1,mo2,mo3,mo3), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 21) >>> ao2mo.outcore.general(mol, (mo1,mo2,mo3,mo3), 'oh2.h5', compact=False) >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 36) >>> ao2mo.outcore.general(mol, (mo1,mo1,mo2,mo2), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (55, 36) >>> ao2mo.outcore.general(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', dataname='new') >>> view('oh2.h5', 'new') dataset ['eri_mo', 'new'], shape (55, 55) >>> ao2mo.outcore.general(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', intor='int2e_ip1_sph', aosym='s1', comp=3) >>> view('oh2.h5') dataset ['eri_mo', 'new'], shape (3, 100, 100) >>> ao2mo.outcore.general(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> view('oh2.h5') dataset ['eri_mo', 'new'], shape (3, 100, 55)
-
pyscf.ao2mo.outcore.
general_iofree
(mol, mo_coeffs, intor='int2e', aosym='s4', comp=None, max_memory=2000, ioblk_size=256, verbose=2, compact=True)¶ For the given four sets of orbitals, transfer arbitrary spherical AO integrals to MO integrals on the fly. This function is a wrap for
ao2mo.outcore.general()
. It’s not really IO free. The returned MO integrals are held in memory. For backward compatibility, it is used to replace the non-existed function direct.general_iofree.- Args:
- mol
Mole
object AO integrals will be generated in terms of mol._atm, mol._bas, mol._env
- mo_coeffs4-item list of ndarray
Four sets of orbital coefficients, corresponding to the four indices of (ij|kl)
- mol
- Kwargs
- intorstr
Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Components of the integrals, e.g. int2e_ip_sph has 3 components.
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
2D/3D MO-integral array. They may or may not have the permutation symmetry, depending on the given orbitals, and the kwargs compact. If the four sets of orbitals are identical, the MO integrals will at most have 4-fold symmetry.
Examples:
>>> from pyscf import gto >>> from pyscf import ao2mo >>> import h5py >>> def view(h5file, dataname='eri_mo'): ... f5 = h5py.File(h5file, 'r') ... print('dataset %s, shape %s' % (str(f5.keys()), str(f5[dataname].shape))) ... f5.close() >>> mol = gto.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)) >>> mo3 = numpy.random.random((mol.nao_nr(), 6)) >>> mo4 = numpy.random.random((mol.nao_nr(), 4)) >>> eri1 = ao2mo.outcore.general_iofree(mol, (mo1,mo2,mo3,mo4)) >>> print(eri1.shape) (80, 24) >>> eri1 = ao2mo.outcore.general_iofree(mol, (mo1,mo2,mo3,mo3)) >>> print(eri1.shape) (80, 21) >>> eri1 = ao2mo.outcore.general_iofree(mol, (mo1,mo2,mo3,mo3), compact=False) >>> print(eri1.shape) (80, 36) >>> eri1 = ao2mo.outcore.general_iofree(mol, (mo1,mo1,mo1,mo1), intor='int2e_ip1_sph', aosym='s1', comp=3) >>> print(eri1.shape) (3, 100, 100) >>> eri1 = ao2mo.outcore.general_iofree(mol, (mo1,mo1,mo1,mo1), intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> print(eri1.shape) (3, 100, 55)
-
pyscf.ao2mo.outcore.
guess_e1bufsize
(max_memory, ioblk_size, nij_pair, nao_pair, comp)¶
-
pyscf.ao2mo.outcore.
guess_e2bufsize
(ioblk_size, nrows, ncols)¶
-
pyscf.ao2mo.outcore.
guess_shell_ranges
(mol, aosym, max_iobuf, max_aobuf=None, ao_loc=None, compress_diag=True)¶
-
pyscf.ao2mo.outcore.
half_e1
(mol, mo_coeffs, swapfile, intor='int2e', aosym='s4', comp=1, max_memory=2000, ioblk_size=256, verbose=2, compact=True, ao2mopt=None)¶ Half transform arbitrary spherical AO integrals to MO integrals for the given two sets of orbitals
- Args:
- mol
Mole
object AO integrals will be generated in terms of mol._atm, mol._bas, mol._env
- mo_coeffndarray
Transform (ij|kl) with the same set of orbitals.
- swapfilestr or h5py File or h5py Group object
To store the transformed integrals, in HDF5 format. The transformed integrals are saved in blocks.
- mol
- Kwargs
- intorstr
Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Components of the integrals, e.g. int2e_ip_sph has 3 components.
- verboseint
Print level
- max_memoryfloat or int
The maximum size of cache to use (in MB), large cache may not improve performance.
- ioblk_sizefloat or int
The block size for IO, large block size may not improve performance
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- ao2mopt
AO2MOpt
object Precomputed data to improve perfomance
- Returns:
None
-
pyscf.ao2mo.outcore.
iden_coeffs
(mo1, mo2)¶
-
pyscf.ao2mo.outcore.
prange
(start, end, step)¶
pyscf.ao2mo.r_outcore module¶
-
pyscf.ao2mo.r_outcore.
full
(mol, mo_coeff, erifile, dataname='eri_mo', intor='int2e_spinor', aosym='s4', comp=None, max_memory=4000, ioblk_size=256, verbose=2)¶
-
pyscf.ao2mo.r_outcore.
full_iofree
(mol, mo_coeff, intor='int2e_spinor', aosym='s4', comp=None, verbose=2, **kwargs)¶
-
pyscf.ao2mo.r_outcore.
general
(mol, mo_coeffs, erifile, dataname='eri_mo', intor='int2e_spinor', aosym='s4', comp=None, max_memory=4000, ioblk_size=256, verbose=2)¶
-
pyscf.ao2mo.r_outcore.
general_iofree
(mol, mo_coeffs, intor='int2e_spinor', aosym='s4', comp=None, verbose=2, **kwargs)¶
-
pyscf.ao2mo.r_outcore.
guess_e1bufsize
(max_memory, ioblk_size, nij_pair, nao_pair, comp)¶
-
pyscf.ao2mo.r_outcore.
guess_e2bufsize
(ioblk_size, nrows, ncols)¶
-
pyscf.ao2mo.r_outcore.
half_e1
(mol, mo_coeffs, swapfile, intor='int2e_spinor', aosym='s4', comp=None, max_memory=4000, ioblk_size=256, verbose=2, ao2mopt=None)¶
-
pyscf.ao2mo.r_outcore.
iden_coeffs
(mo1, mo2)¶
-
pyscf.ao2mo.r_outcore.
prange
(start, end, step)¶
pyscf.ao2mo.semi_incore module¶
Created on Thu May 17 11:05:22 2018
@author: Bryan Lau
A module that will do on-disk transformation of two electron integrals, and also return specific slices of (o)ccupied and (v)irtual ones needed for post HF
Comparing to the full in-memory transformation (see incore.py) which holds all intermediates in memory, this version uses less memory but performs slow due to IO overhead.
-
pyscf.ao2mo.semi_incore.
general
(eri, mo_coeffs, erifile, dataname='eri_mo', ioblk_size=128, compact=True, verbose=3)¶ For the given four sets of orbitals, transfer arbitrary spherical AO integrals to MO integrals on disk. Args:
eri : 8-fold reduced eri vector mo_coeffs : 4-item list of ndarray
Four sets of orbital coefficients, corresponding to the four indices of (ij|kl)
- erifilestr or h5py File or h5py Group object
To store the transformed integrals, in HDF5 format.
- Kwargs
- datanamestr
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 dataname, the new integrals data will overwrite the old one.
- ioblk_sizefloat or int
The block size for IO, large block size may not improve performance
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Pseudocode / algorithm:
u = mu v = nu l = lambda o = sigma
Assume eri’s are 8-fold reduced. nij/nkl_pair = npair or i*j/k*l if only transforming a subset
- First half transform:
- Initialize half_eri of size (nij_pair,npair)
- For lo = 1 -> npair
Unpack row lo Unpack row lo to matrix E_{uv}^{lo} Transform C_ui^+*E*C_nj -> E_{ij}^{lo} Ravel or pack E_{ij}^{lo} Save E_{ij}^{lo} -> half_eri[:,lo]
- Second half transform:
- Initialize h5d_eri of size (nij_pair,nkl_pair)
- For ij = 1 -> nij_pair
Load and unpack half_eri[ij,:] -> E_{lo}^{ij} Transform C_{lk}E_{lo}^{ij}C_{ol} -> E_{kl}^{ij} Repack E_{kl}^{ij} Save E_{kl}^{ij} -> h5d_eri[ij,:]
Each matrix is indexed by the composite index ij x kl, where ij/kl is either npair or ixj/kxl, if only a subset of MOs are being transformed. Since entire rows or columns need to be read in, the arrays are chunked such that IOBLK_SIZE = row/col x chunking col/row. For example, for the first half transform, we would save in nij_pair x IOBLK_SIZE/nij_pair, then load in IOBLK_SIZE/nkl_pair x npair for the second half transform.
—— kl —–> |jxl | ij | | v
As a first guess, the chunking size is jxl. If the super-rows/cols are larger than IOBLK_SIZE, then the chunk rectangle jxl is trimmed accordingly. The pathological limiting case is where the dimensions nao_pair, nij_pair, or nkl_pair are so large that the arrays are chunked 1x1, in which case IOBLK_SIZE needs to be increased.
Module contents¶
General Integral transformation module¶
Simple usage:
>>> from pyscf import gto, scf, ao2mo
>>> mol = gto.M(atom='H 0 0 0; F 0 0 1')
>>> mf = scf.RHF(mol).run()
>>> mo_ints = ao2mo.kernel(mol, mf.mo_coeff)
-
pyscf.ao2mo.
full
(eri_or_mol, mo_coeff, erifile=None, dataname='eri_mo', intor='int2e', *args, **kwargs)¶ MO integral transformation. The four indices (ij|kl) are transformed with the same set of orbitals.
- Args:
- eri_or_molndarray or Mole object
If AO integrals are given as ndarray, it can be either 8-fold or 4-fold symmetry. The integral transformation are computed incore (ie all intermediate are held in memory). If Mole object is given, AO integrals are generated on the fly and outcore algorithm is used (ie intermediate data are held on disk).
- mo_coeffndarray
Orbital coefficients in 2D array
- erifilestr or h5py File or h5py Group object
Note this argument is effective when eri_or_mol is Mole object. The file to store the transformed integrals. If not given, the transformed integrals are held in memory.
- Kwargs:
- datanamestr
Note this argument is effective when eri_or_mol is Mole object. 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 dataname, the new integrals data will overwrite the old one.
- intorstr
Note this argument is effective when eri_or_mol is Mole object. Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Note this argument is effective when eri_or_mol is Mole object. Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Note this argument is effective when eri_or_mol is Mole object. Components of the integrals, e.g. int2e_ip_sph has 3 components.
- max_memoryfloat or int
Note this argument is effective when eri_or_mol is Mole object. The maximum size of cache to use (in MB), large cache may not improve performance.
- ioblk_sizefloat or int
Note this argument is effective when eri_or_mol is Mole object. The block size for IO, large block size may not improve performance
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
If eri_or_mol is array or erifile is not give, the function returns 2D array (or 3D array if comp > 1) of transformed MO integrals. The MO integrals may or may not have the permutation symmetry (controlled by the kwargs compact). Otherwise, return the file/fileobject where the MO integrals are saved.
Examples:
>>> from pyscf import gto, ao2mo >>> import h5py >>> def view(h5file, dataname='eri_mo'): ... with h5py.File(h5file, 'r') as f5: ... print('dataset %s, shape %s' % (str(f5.keys()), str(f5[dataname].shape))) >>> mol = gto.M(atom='O 0 0 0; H 0 1 0; H 0 0 1', basis='sto3g') >>> mo1 = numpy.random.random((mol.nao_nr(), 10))
>>> eri1 = ao2mo.full(mol, mo1) >>> print(eri1.shape) (55, 55)
>>> eri = mol.intor('int2e_sph', aosym='s8') >>> eri1 = ao2mo.full(eri, mo1, compact=False) >>> print(eri1.shape) (100, 100)
>>> ao2mo.full(mol, mo1, 'full.h5') >>> view('full.h5') dataset ['eri_mo'], shape (55, 55)
>>> ao2mo.full(mol, mo1, 'full.h5', dataname='new', compact=False) >>> view('full.h5', 'new') dataset ['eri_mo', 'new'], shape (100, 100)
>>> ao2mo.full(mol, mo1, 'full.h5', intor='int2e_ip1_sph', aosym='s1', comp=3) >>> view('full.h5') dataset ['eri_mo', 'new'], shape (3, 100, 100)
>>> ao2mo.full(mol, mo1, 'full.h5', intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> view('full.h5') dataset ['eri_mo', 'new'], shape (3, 100, 55)
-
pyscf.ao2mo.
general
(eri_or_mol, mo_coeffs, erifile=None, dataname='eri_mo', intor='int2e', *args, **kwargs)¶ Given four sets of orbitals corresponding to the four MO indices, transfer arbitrary spherical AO integrals to MO integrals.
- Args:
- eri_or_molndarray or Mole object
If AO integrals are given as ndarray, it can be either 8-fold or 4-fold symmetry. The integral transformation are computed incore (ie all intermediate are held in memory). If Mole object is given, AO integrals are generated on the fly and outcore algorithm is used (ie intermediate data are held on disk).
- mo_coeffs4-item list of ndarray
Four sets of orbital coefficients, corresponding to the four indices of (ij|kl)
- erifilestr or h5py File or h5py Group object
Note this argument is effective when eri_or_mol is Mole object. The file to store the transformed integrals. If not given, the transformed integrals are held in memory.
- Kwargs:
- datanamestr
Note this argument is effective when eri_or_mol is Mole object. 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 dataname, the new integrals data will overwrite the old one.
- intorstr
Note this argument is effective when eri_or_mol is Mole object. Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Note this argument is effective when eri_or_mol is Mole object. Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Note this argument is effective when eri_or_mol is Mole object. Components of the integrals, e.g. int2e_ip_sph has 3 components.
- max_memoryfloat or int
Note this argument is effective when eri_or_mol is Mole object. The maximum size of cache to use (in MB), large cache may not improve performance.
- ioblk_sizefloat or int
Note this argument is effective when eri_or_mol is Mole object. The block size for IO, large block size may not improve performance
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
If eri_or_mol is array or erifile is not give, the function returns 2D array (or 3D array, if comp > 1) of transformed MO integrals. The MO integrals may at most have 4-fold symmetry (if the four sets of orbitals are identical) or may not have the permutation symmetry (controlled by the kwargs compact). Otherwise, return the file/fileobject where the MO integrals are saved.
Examples:
>>> from pyscf import gto, ao2mo >>> import h5py >>> def view(h5file, dataname='eri_mo'): ... with h5py.File(h5file, 'r') as f5: ... print('dataset %s, shape %s' % (str(f5.keys()), str(f5[dataname].shape))) >>> mol = gto.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)) >>> mo3 = numpy.random.random((mol.nao_nr(), 6)) >>> mo4 = numpy.random.random((mol.nao_nr(), 4))
>>> eri1 = ao2mo.general(eri, (mo1,mo2,mo3,mo4)) >>> print(eri1.shape) (80, 24)
>>> eri1 = ao2mo.general(eri, (mo1,mo2,mo3,mo3)) >>> print(eri1.shape) (80, 21)
>>> eri1 = ao2mo.general(eri, (mo1,mo2,mo3,mo3), compact=False) >>> print(eri1.shape) (80, 36)
>>> eri1 = ao2mo.general(eri, (mo1,mo1,mo2,mo2)) >>> print(eri1.shape) (55, 36)
>>> eri1 = ao2mo.general(eri, (mo1,mo2,mo1,mo2)) >>> print(eri1.shape) (80, 80)
>>> ao2mo.general(mol, (mo1,mo2,mo3,mo4), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 24)
>>> ao2mo.general(mol, (mo1,mo2,mo3,mo3), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 21)
>>> ao2mo.general(mol, (mo1,mo2,mo3,mo3), 'oh2.h5', compact=False) >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 36)
>>> ao2mo.general(mol, (mo1,mo1,mo2,mo2), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (55, 36)
>>> ao2mo.general(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', dataname='new') >>> view('oh2.h5', 'new') dataset ['eri_mo', 'new'], shape (55, 55)
>>> ao2mo.general(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', intor='int2e_ip1_sph', aosym='s1', comp=3) >>> view('oh2.h5') dataset ['eri_mo', 'new'], shape (3, 100, 100)
>>> ao2mo.general(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> view('oh2.h5') dataset ['eri_mo', 'new'], shape (3, 100, 55)
-
pyscf.ao2mo.
get_ao_eri
(mol)¶ 2-electron integrals in AO basis
-
pyscf.ao2mo.
get_mo_eri
(eri_or_mol, mo_coeffs, erifile=None, dataname='eri_mo', intor='int2e', *args, **kwargs)¶ Transfer arbitrary spherical AO integrals to MO integrals, for given orbitals or four sets of orbitals. See also
ao2mo.full()
andao2mo.general()
.- Args:
- eri_or_molndarray or Mole object
If AO integrals are given as ndarray, it can be either 8-fold or 4-fold symmetry. The integral transformation are computed incore (ie all intermediate are held in memory). If Mole object is given, AO integrals are generated on the fly and outcore algorithm is used (ie intermediate data are held on disk).
- mo_coeffsan np array or a list of arrays
A matrix of orbital coefficients if it is a numpy ndarray; Or four sets of orbital coefficients if it is a list of arrays, corresponding to the four indices of (ij|kl).
- Kwargs:
- erifilestr or h5py File or h5py Group object
Note this argument is effective when eri_or_mol is Mole object. The file 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 when eri_or_mol is Mole object. 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 dataname, the new integrals data will overwrite the old one.
- intorstr
Note this argument is effective when eri_or_mol is Mole object. Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Note this argument is effective when eri_or_mol is Mole object. Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Note this argument is effective when eri_or_mol is Mole object. Components of the integrals, e.g. int2e_ip_sph has 3 components.
- max_memoryfloat or int
Note this argument is effective when eri_or_mol is Mole object. The maximum size of cache to use (in MB), large cache may not improve performance.
- ioblk_sizefloat or int
Note this argument is effective when eri_or_mol is Mole object. The block size for IO, large block size may not improve performance
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
If eri_or_mol is array or erifile is not give, the function returns 2D array (or 3D array, if comp > 1) of transformed MO integrals. The MO integrals may at most have 4-fold symmetry (if the four sets of orbitals are identical) or may not have the permutation symmetry (controlled by the kwargs compact). Otherwise, return the file/fileobject where the MO integrals are saved.
Examples:
>>> from pyscf import gto, ao2mo >>> import h5py >>> def view(h5file, dataname='eri_mo'): ... with h5py.File(h5file) as f5: ... print('dataset %s, shape %s' % (str(f5.keys()), str(f5[dataname].shape))) >>> mol = gto.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)) >>> mo3 = numpy.random.random((mol.nao_nr(), 6)) >>> mo4 = numpy.random.random((mol.nao_nr(), 4))
>>> eri1 = ao2mo.kernel(mol, mo1) >>> print(eri1.shape) (55, 55)
>>> eri = mol.intor('int2e_sph', aosym='s8') >>> eri1 = ao2mo.kernel(eri, mo1, compact=False) >>> print(eri1.shape) (100, 100)
>>> ao2mo.kernel(mol, mo1, erifile='full.h5') >>> view('full.h5') dataset ['eri_mo'], shape (55, 55)
>>> ao2mo.kernel(mol, mo1, 'full.h5', dataname='new', compact=False) >>> view('full.h5', 'new') dataset ['eri_mo', 'new'], shape (100, 100)
>>> ao2mo.kernel(mol, mo1, 'full.h5', intor='int2e_ip1_sph', aosym='s1', comp=3) >>> view('full.h5') dataset ['eri_mo', 'new'], shape (3, 100, 100)
>>> ao2mo.kernel(mol, mo1, 'full.h5', intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> view('full.h5') dataset ['eri_mo', 'new'], shape (3, 100, 55)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo2,mo3,mo4)) >>> print(eri1.shape) (80, 24)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo2,mo3,mo3)) >>> print(eri1.shape) (80, 21)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo2,mo3,mo3), compact=False) >>> print(eri1.shape) (80, 36)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo1,mo2,mo2)) >>> print(eri1.shape) (55, 36)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo2,mo1,mo2)) >>> print(eri1.shape) (80, 80)
>>> ao2mo.kernel(mol, (mo1,mo2,mo3,mo4), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 24)
>>> ao2mo.kernel(mol, (mo1,mo2,mo3,mo3), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 21)
>>> ao2mo.kernel(mol, (mo1,mo2,mo3,mo3), 'oh2.h5', compact=False) >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 36)
>>> ao2mo.kernel(mol, (mo1,mo1,mo2,mo2), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (55, 36)
>>> ao2mo.kernel(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', dataname='new') >>> view('oh2.h5', 'new') dataset ['eri_mo', 'new'], shape (55, 55)
>>> ao2mo.kernel(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', intor='int2e_ip1_sph', aosym='s1', comp=3) >>> view('oh2.h5') dataset ['eri_mo', 'new'], shape (3, 100, 100)
>>> ao2mo.kernel(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> view('oh2.h5') dataset ['eri_mo', 'new'], shape (3, 100, 55)
-
pyscf.ao2mo.
kernel
(eri_or_mol, mo_coeffs, erifile=None, dataname='eri_mo', intor='int2e', *args, **kwargs)¶ Transfer arbitrary spherical AO integrals to MO integrals, for given orbitals or four sets of orbitals. See also
ao2mo.full()
andao2mo.general()
.- Args:
- eri_or_molndarray or Mole object
If AO integrals are given as ndarray, it can be either 8-fold or 4-fold symmetry. The integral transformation are computed incore (ie all intermediate are held in memory). If Mole object is given, AO integrals are generated on the fly and outcore algorithm is used (ie intermediate data are held on disk).
- mo_coeffsan np array or a list of arrays
A matrix of orbital coefficients if it is a numpy ndarray; Or four sets of orbital coefficients if it is a list of arrays, corresponding to the four indices of (ij|kl).
- Kwargs:
- erifilestr or h5py File or h5py Group object
Note this argument is effective when eri_or_mol is Mole object. The file 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 when eri_or_mol is Mole object. 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 dataname, the new integrals data will overwrite the old one.
- intorstr
Note this argument is effective when eri_or_mol is Mole object. Name of the 2-electron integral. Ref to
getints_by_shell()
for the complete list of available 2-electron integral names- aosymint or str
Note this argument is effective when eri_or_mol is Mole object. Permutation symmetry for the AO 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‘a4ij’ : 4-fold symmetry with anti-symmetry between i, j in (ij|kl) (TODO)‘a4kl’ : 4-fold symmetry with anti-symmetry between k, l in (ij|kl) (TODO)‘a2ij’ : anti-symmetry between i, j in (ij|kl) (TODO)‘a2kl’ : anti-symmetry between k, l in (ij|kl) (TODO)- compint
Note this argument is effective when eri_or_mol is Mole object. Components of the integrals, e.g. int2e_ip_sph has 3 components.
- max_memoryfloat or int
Note this argument is effective when eri_or_mol is Mole object. The maximum size of cache to use (in MB), large cache may not improve performance.
- ioblk_sizefloat or int
Note this argument is effective when eri_or_mol is Mole object. The block size for IO, large block size may not improve performance
- verboseint
Print level
- compactbool
When compact is True, depending on the four oribital sets, the returned MO integrals has (up to 4-fold) permutation symmetry. If it’s False, the function will abandon any permutation symmetry, and return the “plain” MO integrals
- Returns:
If eri_or_mol is array or erifile is not give, the function returns 2D array (or 3D array, if comp > 1) of transformed MO integrals. The MO integrals may at most have 4-fold symmetry (if the four sets of orbitals are identical) or may not have the permutation symmetry (controlled by the kwargs compact). Otherwise, return the file/fileobject where the MO integrals are saved.
Examples:
>>> from pyscf import gto, ao2mo >>> import h5py >>> def view(h5file, dataname='eri_mo'): ... with h5py.File(h5file) as f5: ... print('dataset %s, shape %s' % (str(f5.keys()), str(f5[dataname].shape))) >>> mol = gto.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)) >>> mo3 = numpy.random.random((mol.nao_nr(), 6)) >>> mo4 = numpy.random.random((mol.nao_nr(), 4))
>>> eri1 = ao2mo.kernel(mol, mo1) >>> print(eri1.shape) (55, 55)
>>> eri = mol.intor('int2e_sph', aosym='s8') >>> eri1 = ao2mo.kernel(eri, mo1, compact=False) >>> print(eri1.shape) (100, 100)
>>> ao2mo.kernel(mol, mo1, erifile='full.h5') >>> view('full.h5') dataset ['eri_mo'], shape (55, 55)
>>> ao2mo.kernel(mol, mo1, 'full.h5', dataname='new', compact=False) >>> view('full.h5', 'new') dataset ['eri_mo', 'new'], shape (100, 100)
>>> ao2mo.kernel(mol, mo1, 'full.h5', intor='int2e_ip1_sph', aosym='s1', comp=3) >>> view('full.h5') dataset ['eri_mo', 'new'], shape (3, 100, 100)
>>> ao2mo.kernel(mol, mo1, 'full.h5', intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> view('full.h5') dataset ['eri_mo', 'new'], shape (3, 100, 55)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo2,mo3,mo4)) >>> print(eri1.shape) (80, 24)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo2,mo3,mo3)) >>> print(eri1.shape) (80, 21)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo2,mo3,mo3), compact=False) >>> print(eri1.shape) (80, 36)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo1,mo2,mo2)) >>> print(eri1.shape) (55, 36)
>>> eri1 = ao2mo.kernel(eri, (mo1,mo2,mo1,mo2)) >>> print(eri1.shape) (80, 80)
>>> ao2mo.kernel(mol, (mo1,mo2,mo3,mo4), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 24)
>>> ao2mo.kernel(mol, (mo1,mo2,mo3,mo3), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 21)
>>> ao2mo.kernel(mol, (mo1,mo2,mo3,mo3), 'oh2.h5', compact=False) >>> view('oh2.h5') dataset ['eri_mo'], shape (80, 36)
>>> ao2mo.kernel(mol, (mo1,mo1,mo2,mo2), 'oh2.h5') >>> view('oh2.h5') dataset ['eri_mo'], shape (55, 36)
>>> ao2mo.kernel(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', dataname='new') >>> view('oh2.h5', 'new') dataset ['eri_mo', 'new'], shape (55, 55)
>>> ao2mo.kernel(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', intor='int2e_ip1_sph', aosym='s1', comp=3) >>> view('oh2.h5') dataset ['eri_mo', 'new'], shape (3, 100, 100)
>>> ao2mo.kernel(mol, (mo1,mo1,mo1,mo1), 'oh2.h5', intor='int2e_ip1_sph', aosym='s2kl', comp=3) >>> view('oh2.h5') dataset ['eri_mo', 'new'], shape (3, 100, 55)