7. ao2mo — Integral transformations

7.1. 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)

7.2. incore

pyscf.ao2mo.incore.full(eri_ao, mo_coeff, verbose=0, compact=True, **kwargs)[source]

MO integral transformation for the given orbital.

Args:
eri_ao
: ndarray
AO integrals, can be either 8-fold or 4-fold symmetry.
mo_coeff
: ndarray
Transform (ij|kl) with the same set of orbitals.
Kwargs:
verbose
: int
Print level
compact
: bool
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)[source]

For the given four sets of orbitals, transfer the 8-fold or 4-fold 2e AO integrals to MO integrals.

Args:
eri_ao
: ndarray
AO integrals, can be either 8-fold or 4-fold symmetry.
mo_coeffs
: 4-item list of ndarray
Four sets of orbital coefficients, corresponding to the four indices of (ij|kl)
Kwargs:
verbose
: int
Print level
compact
: bool
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)[source]

Given two set of orbitals, half transform the (ij| pair of 8-fold or 4-fold AO integrals (ij|kl)

Args:
eri_ao
: ndarray
AO integrals, can be either 8-fold or 4-fold symmetry.
mo_coeffs
: list of ndarray
Two sets of orbital coefficients, corresponding to the i, j indices of (ij|kl)
Kwargs:
compact
: bool
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)

7.3. outcore

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

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_coeff
: ndarray
Transform (ij|kl) with the same set of orbitals.
erifile
: str or h5py File or h5py Group object
To store the transformed integrals, in HDF5 format.
Kwargs:
dataname
: str
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.
intor
: str
Name of the 2-electron integral. Ref to getints_by_shell() for the complete list of available 2-electron integral names
aosym
: int 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)
comp
: int
Components of the integrals, e.g. int2e_ip_sph has 3 components.
max_memory
: float or int
The maximum size of cache to use (in MB), large cache may not improve performance.
ioblk_size
: float or int
The block size for IO, large block size may not improve performance
verbose
: int
Print level
compact
: bool
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)
...     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)[source]

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_coeff
: ndarray
Transform (ij|kl) with the same set of orbitals.
erifile
: str
To store the transformed integrals, in HDF5 format.
Kwargs
dataname
: str
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.
intor
: str
Name of the 2-electron integral. Ref to getints_by_shell() for the complete list of available 2-electron integral names
aosym
: int 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)
comp
: int
Components of the integrals, e.g. int2e_ip_sph has 3 components.
verbose
: int
Print level
max_memory
: float or int
The maximum size of cache to use (in MB), large cache may not improve performance.
ioblk_size
: float or int
The block size for IO, large block size may not improve performance
verbose
: int
Print level
compact
: bool
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)[source]

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_coeffs
: 4-item list of ndarray
Four sets of orbital coefficients, corresponding to the four indices of (ij|kl)
erifile
: str or h5py File or h5py Group object
To store the transformed integrals, in HDF5 format.
Kwargs
dataname
: str
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.
intor
: str
Name of the 2-electron integral. Ref to getints_by_shell() for the complete list of available 2-electron integral names
aosym
: int 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)
comp
: int
Components of the integrals, e.g. int2e_ip_sph has 3 components.
max_memory
: float or int
The maximum size of cache to use (in MB), large cache may not improve performance.
ioblk_size
: float or int
The block size for IO, large block size may not improve performance
verbose
: int
Print level
compact
: bool
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)
...     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)[source]

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_coeffs
: 4-item list of ndarray
Four sets of orbital coefficients, corresponding to the four indices of (ij|kl)
Kwargs
intor
: str
Name of the 2-electron integral. Ref to getints_by_shell() for the complete list of available 2-electron integral names
aosym
: int 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)
comp
: int
Components of the integrals, e.g. int2e_ip_sph has 3 components.
verbose
: int
Print level
compact
: bool
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)
...     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.half_e1(mol, mo_coeffs, swapfile, intor='int2e', aosym='s4', comp=1, max_memory=2000, ioblk_size=256, verbose=2, compact=True, ao2mopt=None)[source]

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_coeff
: ndarray
Transform (ij|kl) with the same set of orbitals.
swapfile
: str or h5py File or h5py Group object
To store the transformed integrals, in HDF5 format. The transformed integrals are saved in blocks.
Kwargs
intor
: str
Name of the 2-electron integral. Ref to getints_by_shell() for the complete list of available 2-electron integral names
aosym
: int 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)
comp
: int
Components of the integrals, e.g. int2e_ip_sph has 3 components.
verbose
: int
Print level
max_memory
: float or int
The maximum size of cache to use (in MB), large cache may not improve performance.
ioblk_size
: float or int
The block size for IO, large block size may not improve performance
verbose
: int
Print level
compact
: bool
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

7.4. addons

class pyscf.ao2mo.addons.load(eri, dataname='eri_mo')[source]

load 2e integrals from hdf5 file

Usage:
with load(erifile) as eri:
print(eri.shape)
pyscf.ao2mo.addons.restore(symmetry, eri, norb, tao=None)[source]

Convert the 2e integrals (in Chemist’s notation) between different level of permutation symmetry (8-fold, 4-fold, or no symmetry)

Args:
symmetry
: int 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.

eri
: ndarray
The symmetry of eri is determined by the size of eri and norb
norb
: int
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)/2
4 : (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,)