pyscf.lib package¶
Submodules¶
pyscf.lib.chkfile module¶
-
pyscf.lib.chkfile.
dump
(chkfile, key, value)¶ Save array(s) in chkfile
- Args:
- chkfilestr
Name of chkfile.
- keystr
key to be used in h5py object. It can contain “/” to represent the path in the HDF5 storage structure.
- valuearray, vector, list … or dict
If value is a python dict or list, the key/value of the dict will be saved recursively as the HDF5 group/dataset structure.
- Returns:
No return value
Examples:
>>> import h5py >>> from pyscf import lib >>> ci = {'Ci' : {'op': ('E', 'i'), 'irrep': ('Ag', 'Au')}} >>> lib.chkfile.save('symm.chk', 'symm', ci) >>> f = h5py.File('symm.chk', 'r') >>> f.keys() ['symm'] >>> f['symm'].keys() ['Ci'] >>> f['symm/Ci'].keys() ['op', 'irrep'] >>> f['symm/Ci/op'] <HDF5 dataset "op": shape (2,), type "|S1">
-
pyscf.lib.chkfile.
dump_chkfile_key
(chkfile, key, value)¶ Save array(s) in chkfile
- Args:
- chkfilestr
Name of chkfile.
- keystr
key to be used in h5py object. It can contain “/” to represent the path in the HDF5 storage structure.
- valuearray, vector, list … or dict
If value is a python dict or list, the key/value of the dict will be saved recursively as the HDF5 group/dataset structure.
- Returns:
No return value
Examples:
>>> import h5py >>> from pyscf import lib >>> ci = {'Ci' : {'op': ('E', 'i'), 'irrep': ('Ag', 'Au')}} >>> lib.chkfile.save('symm.chk', 'symm', ci) >>> f = h5py.File('symm.chk', 'r') >>> f.keys() ['symm'] >>> f['symm'].keys() ['Ci'] >>> f['symm/Ci'].keys() ['op', 'irrep'] >>> f['symm/Ci/op'] <HDF5 dataset "op": shape (2,), type "|S1">
-
pyscf.lib.chkfile.
dump_mol
(mol, chkfile)¶ Save Mole object in chkfile
- Args:
- chkfile str:
Name of chkfile.
- Returns:
No return value
-
pyscf.lib.chkfile.
load
(chkfile, key)¶ Load array(s) from chkfile
- Args:
- chkfilestr
Name of chkfile. The chkfile needs to be saved in HDF5 format.
- keystr
HDF5.dataset name or group name. If key is the name of a HDF5 group, the group will be loaded into a Python dict, recursively.
- Returns:
whatever read from chkfile
Examples:
>>> from pyscf import gto, scf, lib >>> mol = gto.M(atom='He 0 0 0') >>> mf = scf.RHF(mol) >>> mf.chkfile = 'He.chk' >>> mf.kernel() >>> mo_coeff = lib.chkfile.load('He.chk', 'scf/mo_coeff') >>> mo_coeff.shape (1, 1) >>> scfdat = lib.chkfile.load('He.chk', 'scf') >>> scfdat.keys() ['e_tot', 'mo_occ', 'mo_energy', 'mo_coeff']
-
pyscf.lib.chkfile.
load_chkfile_key
(chkfile, key)¶ Load array(s) from chkfile
- Args:
- chkfilestr
Name of chkfile. The chkfile needs to be saved in HDF5 format.
- keystr
HDF5.dataset name or group name. If key is the name of a HDF5 group, the group will be loaded into a Python dict, recursively.
- Returns:
whatever read from chkfile
Examples:
>>> from pyscf import gto, scf, lib >>> mol = gto.M(atom='He 0 0 0') >>> mf = scf.RHF(mol) >>> mf.chkfile = 'He.chk' >>> mf.kernel() >>> mo_coeff = lib.chkfile.load('He.chk', 'scf/mo_coeff') >>> mo_coeff.shape (1, 1) >>> scfdat = lib.chkfile.load('He.chk', 'scf') >>> scfdat.keys() ['e_tot', 'mo_occ', 'mo_energy', 'mo_coeff']
-
pyscf.lib.chkfile.
load_mol
(chkfile)¶ Load Mole object from chkfile. The save_mol/load_mol operation can be used a serialization method for Mole object.
- Args:
- chkfilestr
Name of chkfile.
- Returns:
A (initialized/built) Mole object
Examples:
>>> from pyscf import gto, lib >>> mol = gto.M(atom='He 0 0 0') >>> lib.chkfile.save_mol(mol, 'He.chk') >>> lib.chkfile.load_mol('He.chk') <pyscf.gto.mole.Mole object at 0x7fdcd94d7f50>
-
pyscf.lib.chkfile.
save
(chkfile, key, value)¶ Save array(s) in chkfile
- Args:
- chkfilestr
Name of chkfile.
- keystr
key to be used in h5py object. It can contain “/” to represent the path in the HDF5 storage structure.
- valuearray, vector, list … or dict
If value is a python dict or list, the key/value of the dict will be saved recursively as the HDF5 group/dataset structure.
- Returns:
No return value
Examples:
>>> import h5py >>> from pyscf import lib >>> ci = {'Ci' : {'op': ('E', 'i'), 'irrep': ('Ag', 'Au')}} >>> lib.chkfile.save('symm.chk', 'symm', ci) >>> f = h5py.File('symm.chk', 'r') >>> f.keys() ['symm'] >>> f['symm'].keys() ['Ci'] >>> f['symm/Ci'].keys() ['op', 'irrep'] >>> f['symm/Ci/op'] <HDF5 dataset "op": shape (2,), type "|S1">
-
pyscf.lib.chkfile.
save_mol
(mol, chkfile)¶ Save Mole object in chkfile
- Args:
- chkfile str:
Name of chkfile.
- Returns:
No return value
pyscf.lib.diis module¶
DIIS
-
class
pyscf.lib.diis.
DIIS
(dev=None, filename=None, incore=False)¶ Bases:
object
Direct inversion in the iterative subspace method.
- Attributes:
- spaceint
DIIS subspace size. The maximum number of the vectors to be stored.
- min_space
The minimal size of subspace before DIIS extrapolation.
- Functions:
- update(x, xerr=None) :
If xerr the error vector is given, this function will push the target vector and error vector in the DIIS subspace, and use the error vector to extrapolate the vector and return the extrapolated vector. If xerr is None, this function will take the difference between the current given vector and the last given vector as the error vector to extrapolate the vector.
Examples:
>>> from pyscf import gto, scf, lib >>> mol = gto.M(atom='H 0 0 0; H 0 0 1', basis='ccpvdz') >>> mf = scf.RHF(mol) >>> h = mf.get_hcore() >>> s = mf.get_ovlp() >>> e, c = mf.eig(h, s) >>> occ = mf.get_occ(e, c) >>> # DIIS without error vector >>> adiis = lib.diis.DIIS() >>> for i in range(7): ... dm = mf.make_rdm1(c, occ) ... f = h + mf.get_veff(mol, dm) ... if i > 1: ... f = adiis.update(f) ... e, c = mf.eig(f, s) ... print('E_%d = %.12f' % (i, mf.energy_tot(dm, h, mf.get_veff(mol, dm)))) E_0 = -1.050329433306 E_1 = -1.098566175145 E_2 = -1.100103795287 E_3 = -1.100152104615 E_4 = -1.100153706922 E_5 = -1.100153764848 E_6 = -1.100153764878
>>> # Take Hartree-Fock gradients as the error vector >>> adiis = lib.diis.DIIS() >>> for i in range(7): ... dm = mf.make_rdm1(c, occ) ... f = h + mf.get_veff(mol, dm) ... if i > 1: ... f = adiis.update(f, mf.get_grad(c, occ, f)) ... e, c = mf.eig(f, s) ... print('E_%d = %.12f' % (i, mf.energy_tot(dm, h, mf.get_veff(mol, dm)))) E_0 = -1.050329433306 E_1 = -1.098566175145 E_2 = -1.100103795287 E_3 = -1.100152104615 E_4 = -1.100153763813 E_5 = -1.100153764878 E_6 = -1.100153764878
-
extrapolate
(nd=None)¶
-
get_err_vec
(idx)¶
-
get_num_vec
()¶
-
get_vec
(idx)¶
-
push_err_vec
(xerr)¶
-
push_vec
(x)¶
-
restore
(filename, inplace=True)¶ Read diis contents from a diis file and replace the attributes of current diis object if needed, then construct the vector.
-
update
(x, xerr=None)¶ Extrapolate vector
If xerr the error vector is given, this function will push the target
vector and error vector in the DIIS subspace, and use the error vector to extrapolate the vector and return the extrapolated vector. * If xerr is None, this function will take the difference between the current given vector and the last given vector as the error vector to extrapolate the vector.
-
pyscf.lib.diis.
restore
(filename)¶ Restore/construct diis object based on a diis file
pyscf.lib.exceptions module¶
Errors and exceptions
-
exception
pyscf.lib.exceptions.
BasisNotFoundError
¶ Bases:
RuntimeError
pyscf.lib.linalg_helper module¶
Extension to scipy.linalg module
-
pyscf.lib.linalg_helper.
cho_solve
(a, b, strict_sym_pos=True)¶ Solve ax = b, where a is a postive definite hermitian matrix
- Kwargs:
- strict_sym_pos (bool)Whether to impose the strict positive definition
on matrix a
-
pyscf.lib.linalg_helper.
davidson
(aop, x0, precond, tol=1e-12, max_cycle=50, max_space=12, lindep=1e-14, max_memory=2000, dot=<function dot>, callback=None, nroots=1, lessio=False, pick=None, verbose=2, follow_state=False)¶ Davidson diagonalization method to solve a c = e c. Ref [1] E.R. Davidson, J. Comput. Phys. 17 (1), 87-94 (1975). [2] http://people.inf.ethz.ch/arbenz/ewp/Lnotes/chapter11.pdf
Note: This function has an overhead of memory usage ~4*x0.size*nroots
- Args:
- aopfunction(x) => array_like_x
aop(x) to mimic the matrix vector multiplication \(\sum_{j}a_{ij}*x_j\). The argument is a 1D array. The returned value is a 1D array.
- x01D array or a list of 1D array
Initial guess. The initial guess vector(s) are just used as the initial subspace bases. If the subspace is smaller than “nroots”, eg 10 roots and one initial guess, all eigenvectors are chosen as the eigenvectors during the iterations. The first iteration has one eigenvector, the next iteration has two, the third iterastion has 4, …, until the subspace size > nroots.
- preconddiagonal elements of the matrix or function(dx, e, x0) => array_like_dx
Preconditioner to generate new trial vector. The argument dx is a residual vector
a*x0-e*x0
; e is the current eigenvalue; x0 is the current eigenvector.
- Kwargs:
- tolfloat
Convergence tolerance.
- max_cycleint
max number of iterations.
- max_spaceint
space size to hold trial vectors.
- lindepfloat
Linear dependency threshold. The function is terminated when the smallest eigenvalue of the metric of the trial vectors is lower than this threshold.
- max_memoryint or float
Allowed memory in MB.
- dotfunction(x, y) => scalar
Inner product
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current envrionment.- nrootsint
Number of eigenvalues to be computed. When nroots > 1, it affects the shape of the return value
- lessiobool
How to compute a*x0 for current eigenvector x0. There are two ways to compute a*x0. One is to assemble the existed a*x. The other is to call aop(x0). The default is the first method which needs more IO and less computational cost. When IO is slow, the second method can be considered.
- pickfunction(w,v,nroots) => (e[idx], w[:,idx], idx)
Function to filter eigenvalues and eigenvectors.
- follow_statebool
If the solution dramatically changes in two iterations, clean the subspace and restart the iteration with the old solution. It can help to improve numerical stability. Default is False.
- Returns:
- efloat or list of floats
Eigenvalue. By default it’s one float number. If
nroots
> 1, it is a list of floats for the lowestnroots
eigenvalues.- c1D array or list of 1D arrays
Eigenvector. By default it’s a 1D array. If
nroots
> 1, it is a list of arrays for the lowestnroots
eigenvectors.
Examples:
>>> from pyscf import lib >>> a = numpy.random.random((10,10)) >>> a = a + a.T >>> aop = lambda x: numpy.dot(a,x) >>> precond = lambda dx, e, x0: dx/(a.diagonal()-e) >>> x0 = a[0] >>> e, c = lib.davidson(aop, x0, precond)
-
pyscf.lib.linalg_helper.
davidson1
(aop, x0, precond, tol=1e-12, max_cycle=50, max_space=12, lindep=1e-14, max_memory=2000, dot=<function dot>, callback=None, nroots=1, lessio=False, pick=None, verbose=2, follow_state=False, tol_residual=None)¶ Davidson diagonalization method to solve a c = e c. Ref [1] E.R. Davidson, J. Comput. Phys. 17 (1), 87-94 (1975). [2] http://people.inf.ethz.ch/arbenz/ewp/Lnotes/chapter11.pdf
Note: This function has an overhead of memory usage ~4*x0.size*nroots
- Args:
- aopfunction([x]) => [array_like_x]
Matrix vector multiplication \(y_{ki} = \sum_{j}a_{ij}*x_{jk}\).
- x01D array or a list of 1D arrays or a function to generate x0 array(s)
Initial guess. The initial guess vector(s) are just used as the initial subspace bases. If the subspace is smaller than “nroots”, eg 10 roots and one initial guess, all eigenvectors are chosen as the eigenvectors during the iterations. The first iteration has one eigenvector, the next iteration has two, the third iterastion has 4, …, until the subspace size > nroots.
- preconddiagonal elements of the matrix or function(dx, e, x0) => array_like_dx
Preconditioner to generate new trial vector. The argument dx is a residual vector
a*x0-e*x0
; e is the current eigenvalue; x0 is the current eigenvector.
- Kwargs:
- tolfloat
Convergence tolerance.
- max_cycleint
max number of iterations.
- max_spaceint
space size to hold trial vectors.
- lindepfloat
Linear dependency threshold. The function is terminated when the smallest eigenvalue of the metric of the trial vectors is lower than this threshold.
- max_memoryint or float
Allowed memory in MB.
- dotfunction(x, y) => scalar
Inner product
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current envrionment.- nrootsint
Number of eigenvalues to be computed. When nroots > 1, it affects the shape of the return value
- lessiobool
How to compute a*x0 for current eigenvector x0. There are two ways to compute a*x0. One is to assemble the existed a*x. The other is to call aop(x0). The default is the first method which needs more IO and less computational cost. When IO is slow, the second method can be considered.
- pickfunction(w,v,nroots) => (e[idx], w[:,idx], idx)
Function to filter eigenvalues and eigenvectors.
- follow_statebool
If the solution dramatically changes in two iterations, clean the subspace and restart the iteration with the old solution. It can help to improve numerical stability. Default is False.
- Returns:
- convbool
Converged or not
- elist of floats
The lowest
nroots
eigenvalues.- clist of 1D arrays
The lowest
nroots
eigenvectors.
Examples:
>>> from pyscf import lib >>> a = numpy.random.random((10,10)) >>> a = a + a.T >>> aop = lambda xs: [numpy.dot(a,x) for x in xs] >>> precond = lambda dx, e, x0: dx/(a.diagonal()-e) >>> x0 = a[0] >>> e, c = lib.davidson(aop, x0, precond, nroots=2) >>> len(e) 2
-
pyscf.lib.linalg_helper.
davidson_nosym
(aop, x0, precond, tol=1e-12, max_cycle=50, max_space=12, lindep=1e-14, max_memory=2000, dot=<function dot>, callback=None, nroots=1, lessio=False, left=False, pick=<function pick_real_eigs>, verbose=2, follow_state=False)¶ Davidson diagonalization to solve the non-symmetric eigenvalue problem
- Args:
- aopfunction([x]) => [array_like_x]
Matrix vector multiplication \(y_{ki} = \sum_{j}a_{ij}*x_{jk}\).
- x01D array or a list of 1D array
Initial guess. The initial guess vector(s) are just used as the initial subspace bases. If the subspace is smaller than “nroots”, eg 10 roots and one initial guess, all eigenvectors are chosen as the eigenvectors during the iterations. The first iteration has one eigenvector, the next iteration has two, the third iterastion has 4, …, until the subspace size > nroots.
- preconddiagonal elements of the matrix or function(dx, e, x0) => array_like_dx
Preconditioner to generate new trial vector. The argument dx is a residual vector
a*x0-e*x0
; e is the current eigenvalue; x0 is the current eigenvector.
- Kwargs:
- tolfloat
Convergence tolerance.
- max_cycleint
max number of iterations.
- max_spaceint
space size to hold trial vectors.
- lindepfloat
Linear dependency threshold. The function is terminated when the smallest eigenvalue of the metric of the trial vectors is lower than this threshold.
- max_memoryint or float
Allowed memory in MB.
- dotfunction(x, y) => scalar
Inner product
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current envrionment.- nrootsint
Number of eigenvalues to be computed. When nroots > 1, it affects the shape of the return value
- lessiobool
How to compute a*x0 for current eigenvector x0. There are two ways to compute a*x0. One is to assemble the existed a*x. The other is to call aop(x0). The default is the first method which needs more IO and less computational cost. When IO is slow, the second method can be considered.
- leftbool
Whether to calculate and return left eigenvectors. Default is False.
- pickfunction(w,v,nroots) => (e[idx], w[:,idx], idx)
Function to filter eigenvalues and eigenvectors.
- follow_statebool
If the solution dramatically changes in two iterations, clean the subspace and restart the iteration with the old solution. It can help to improve numerical stability. Default is False.
- Returns:
- convbool
Converged or not
- elist of eigenvalues
The eigenvalues can be sorted real or complex, depending on the return value of
pick
function.- vllist of 1D arrays
Left eigenvectors. Only returned if
left=True
.- clist of 1D arrays
Right eigenvectors.
Examples:
>>> from pyscf import lib >>> a = numpy.random.random((10,10)) >>> a = a >>> aop = lambda xs: [numpy.dot(a,x) for x in xs] >>> precond = lambda dx, e, x0: dx/(a.diagonal()-e) >>> x0 = a[0] >>> e, vl, vr = lib.davidson(aop, x0, precond, nroots=2, left=True) >>> len(e) 2
-
pyscf.lib.linalg_helper.
davidson_nosym1
(aop, x0, precond, tol=1e-12, max_cycle=50, max_space=12, lindep=1e-14, max_memory=2000, dot=<function dot>, callback=None, nroots=1, lessio=False, left=False, pick=<function pick_real_eigs>, verbose=2, follow_state=False, tol_residual=None)¶
-
pyscf.lib.linalg_helper.
dgeev
(abop, x0, precond, type=1, tol=1e-12, max_cycle=50, max_space=12, lindep=1e-14, max_memory=2000, dot=<function dot>, callback=None, nroots=1, lessio=False, verbose=2)¶ Davidson diagonalization method to solve A c = e B c.
- Args:
- abopfunction(x) => (array_like_x, array_like_x)
abop applies two matrix vector multiplications and returns tuple (Ax, Bx)
- x01D array
Initial guess
- preconddiagonal elements of the matrix or function(dx, e, x0) => array_like_dx
Preconditioner to generate new trial vector. The argument dx is a residual vector
a*x0-e*x0
; e is the current eigenvalue; x0 is the current eigenvector.
- Kwargs:
- tolfloat
Convergence tolerance.
- max_cycleint
max number of iterations.
- max_spaceint
space size to hold trial vectors.
- lindepfloat
Linear dependency threshold. The function is terminated when the smallest eigenvalue of the metric of the trial vectors is lower than this threshold.
- max_memoryint or float
Allowed memory in MB.
- dotfunction(x, y) => scalar
Inner product
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current envrionment.- nrootsint
Number of eigenvalues to be computed. When nroots > 1, it affects the shape of the return value
- lessiobool
How to compute a*x0 for current eigenvector x0. There are two ways to compute a*x0. One is to assemble the existed a*x. The other is to call aop(x0). The default is the first method which needs more IO and less computational cost. When IO is slow, the second method can be considered.
- Returns:
- elist of floats
The lowest
nroots
eigenvalues.- clist of 1D arrays
The lowest
nroots
eigenvectors.
-
pyscf.lib.linalg_helper.
dgeev1
(abop, x0, precond, type=1, tol=1e-12, max_cycle=50, max_space=12, lindep=1e-14, max_memory=2000, dot=<function dot>, callback=None, nroots=1, lessio=False, verbose=2, tol_residual=None)¶ Davidson diagonalization method to solve A c = e B c.
- Args:
- abopfunction([x]) => ([array_like_x], [array_like_x])
abop applies two matrix vector multiplications and returns tuple (Ax, Bx)
- x01D array
Initial guess
- preconddiagonal elements of the matrix or function(dx, e, x0) => array_like_dx
Preconditioner to generate new trial vector. The argument dx is a residual vector
a*x0-e*x0
; e is the current eigenvalue; x0 is the current eigenvector.
- Kwargs:
- tolfloat
Convergence tolerance.
- max_cycleint
max number of iterations.
- max_spaceint
space size to hold trial vectors.
- lindepfloat
Linear dependency threshold. The function is terminated when the smallest eigenvalue of the metric of the trial vectors is lower than this threshold.
- max_memoryint or float
Allowed memory in MB.
- dotfunction(x, y) => scalar
Inner product
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current envrionment.- nrootsint
Number of eigenvalues to be computed. When nroots > 1, it affects the shape of the return value
- lessiobool
How to compute a*x0 for current eigenvector x0. There are two ways to compute a*x0. One is to assemble the existed a*x. The other is to call aop(x0). The default is the first method which needs more IO and less computational cost. When IO is slow, the second method can be considered.
- Returns:
- convbool
Converged or not
- elist of floats
The lowest
nroots
eigenvalues.- clist of 1D arrays
The lowest
nroots
eigenvectors.
-
pyscf.lib.linalg_helper.
dsolve
(aop, b, precond, tol=1e-12, max_cycle=30, dot=<function dot>, lindep=1e-15, verbose=0, tol_residual=None)¶ Davidson iteration to solve linear equation. It works bad.
-
pyscf.lib.linalg_helper.
dsyev
(a, *args, **kwargs)¶
-
pyscf.lib.linalg_helper.
eig
(aop, x0, precond, tol=1e-12, max_cycle=50, max_space=12, lindep=1e-14, max_memory=2000, dot=<function dot>, callback=None, nroots=1, lessio=False, left=False, pick=<function pick_real_eigs>, verbose=2, follow_state=False)¶ Davidson diagonalization to solve the non-symmetric eigenvalue problem
- Args:
- aopfunction([x]) => [array_like_x]
Matrix vector multiplication \(y_{ki} = \sum_{j}a_{ij}*x_{jk}\).
- x01D array or a list of 1D array
Initial guess. The initial guess vector(s) are just used as the initial subspace bases. If the subspace is smaller than “nroots”, eg 10 roots and one initial guess, all eigenvectors are chosen as the eigenvectors during the iterations. The first iteration has one eigenvector, the next iteration has two, the third iterastion has 4, …, until the subspace size > nroots.
- preconddiagonal elements of the matrix or function(dx, e, x0) => array_like_dx
Preconditioner to generate new trial vector. The argument dx is a residual vector
a*x0-e*x0
; e is the current eigenvalue; x0 is the current eigenvector.
- Kwargs:
- tolfloat
Convergence tolerance.
- max_cycleint
max number of iterations.
- max_spaceint
space size to hold trial vectors.
- lindepfloat
Linear dependency threshold. The function is terminated when the smallest eigenvalue of the metric of the trial vectors is lower than this threshold.
- max_memoryint or float
Allowed memory in MB.
- dotfunction(x, y) => scalar
Inner product
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current envrionment.- nrootsint
Number of eigenvalues to be computed. When nroots > 1, it affects the shape of the return value
- lessiobool
How to compute a*x0 for current eigenvector x0. There are two ways to compute a*x0. One is to assemble the existed a*x. The other is to call aop(x0). The default is the first method which needs more IO and less computational cost. When IO is slow, the second method can be considered.
- leftbool
Whether to calculate and return left eigenvectors. Default is False.
- pickfunction(w,v,nroots) => (e[idx], w[:,idx], idx)
Function to filter eigenvalues and eigenvectors.
- follow_statebool
If the solution dramatically changes in two iterations, clean the subspace and restart the iteration with the old solution. It can help to improve numerical stability. Default is False.
- Returns:
- convbool
Converged or not
- elist of eigenvalues
The eigenvalues can be sorted real or complex, depending on the return value of
pick
function.- vllist of 1D arrays
Left eigenvectors. Only returned if
left=True
.- clist of 1D arrays
Right eigenvectors.
Examples:
>>> from pyscf import lib >>> a = numpy.random.random((10,10)) >>> a = a >>> aop = lambda xs: [numpy.dot(a,x) for x in xs] >>> precond = lambda dx, e, x0: dx/(a.diagonal()-e) >>> x0 = a[0] >>> e, vl, vr = lib.davidson(aop, x0, precond, nroots=2, left=True) >>> len(e) 2
-
pyscf.lib.linalg_helper.
eigh
(a, *args, **kwargs)¶
-
pyscf.lib.linalg_helper.
eigh_by_blocks
(h, s=None, labels=None)¶ Solve an ordinary or generalized eigenvalue problem for diagonal blocks. The diagonal blocks are extracted based on the given basis “labels”. The rows and columns which have the same labels are put in the same block. One common scenario one needs the block-wise diagonalization is to diagonalize the matrix in symmetry adapted basis, in which “labels” is the irreps of each basis.
- Args:
- h, s2D array
Complex Hermitian or real symmetric matrix.
- Kwargs:
labels : list
- Returns:
w, v. w is the eigenvalue vector; v is the eigenfunction array; seig is the eigenvalue vector of the metric s.
Examples:
>>> from pyscf import lib >>> a = numpy.ones((4,4)) >>> a[0::3,0::3] = 0 >>> a[1::3,1::3] = 2 >>> a[2::3,2::3] = 4 >>> labels = ['a', 'b', 'c', 'a'] >>> lib.eigh_by_blocks(a, labels) (array([ 0., 0., 2., 4.]), array([[ 1., 0., 0., 0.], [ 0., 0., 1., 0.], [ 0., 0., 0., 1.], [ 0., 1., 0., 0.]])) >>> numpy.linalg.eigh(a) (array([ -8.82020545e-01, -1.81556477e-16, 1.77653793e+00, 5.10548262e+00]), array([[ 6.40734630e-01, -7.07106781e-01, 1.68598330e-01, -2.47050070e-01], [ -3.80616542e-01, 9.40505244e-17, 8.19944479e-01, -4.27577008e-01], [ -1.84524565e-01, 9.40505244e-17, -5.20423152e-01, -8.33732828e-01], [ 6.40734630e-01, 7.07106781e-01, 1.68598330e-01, -2.47050070e-01]]))
>>> from pyscf import gto, lib, symm >>> mol = gto.M(atom='H 0 0 0; H 0 0 1', basis='ccpvdz', symmetry=True) >>> c = numpy.hstack(mol.symm_orb) >>> vnuc_so = reduce(numpy.dot, (c.T, mol.intor('int1e_nuc_sph'), c)) >>> orbsym = symm.label_orb_symm(mol, mol.irrep_name, mol.symm_orb, c) >>> lib.eigh_by_blocks(vnuc_so, labels=orbsym) (array([-4.50766885, -1.80666351, -1.7808565 , -1.7808565 , -1.74189134, -0.98998583, -0.98998583, -0.40322226, -0.30242374, -0.07608981]), ...)
-
pyscf.lib.linalg_helper.
krylov
(aop, b, x0=None, tol=1e-10, max_cycle=30, dot=<function dot>, lindep=1e-15, callback=None, hermi=False, max_memory=2000, verbose=2)¶ Krylov subspace method to solve (1+a) x = b. Ref: J. A. Pople et al, Int. J. Quantum. Chem. Symp. 13, 225 (1979).
- Args:
- aopfunction(x) => array_like_x
aop(x) to mimic the matrix vector multiplication \(\sum_{j}a_{ij} x_j\). The argument is a 1D array. The returned value is a 1D array.
b : a vector or a list of vectors
- Kwargs:
- x01D array
Initial guess
- tolfloat
Tolerance to terminate the operation aop(x).
- max_cycleint
max number of iterations.
- lindepfloat
Linear dependency threshold. The function is terminated when the smallest eigenvalue of the metric of the trial vectors is lower than this threshold.
- dotfunction(x, y) => scalar
Inner product
- callbackfunction(envs_dict) => None
callback function takes one dict as the argument which is generated by the builtin function
locals()
, so that the callback function can access all local variables in the current envrionment.
- Returns:
x : ndarray like b
Examples:
>>> from pyscf import lib >>> a = numpy.random.random((10,10)) * 1e-2 >>> b = numpy.random.random(10) >>> aop = lambda x: numpy.dot(a,x) >>> x = lib.krylov(aop, b) >>> numpy.allclose(numpy.dot(a,x)+x, b) True
-
pyscf.lib.linalg_helper.
make_diag_precond
(diag, level_shift=0)¶ Generate the preconditioner function with the diagonal function.
-
pyscf.lib.linalg_helper.
pick_real_eigs
(w, v, nroots, envs)¶ This function searchs the real eigenvalues or eigenvalues with small imaginary component.
-
pyscf.lib.linalg_helper.
safe_eigh
(h, s, lindep=1e-15)¶ Solve generalized eigenvalue problem h v = w s v.
Note
The number of eigenvalues and eigenvectors might be less than the matrix dimension if linear dependency is found in metric s.
- Args:
- h, s2D array
Complex Hermitian or real symmetric matrix.
- Kwargs:
- lindepfloat
Linear dependency threshold. By diagonalizing the metric s, we consider the eigenvectors are linearly dependent subsets if their eigenvalues are smaller than this threshold.
- Returns:
w, v, seig. w is the eigenvalue vector; v is the eigenfunction array; seig is the eigenvalue vector of the metric s.
pyscf.lib.logger module¶
Logging system
Log level¶
Level |
number |
DEBUG4 |
9 |
DEBUG3 |
8 |
DEBUG2 |
7 |
DEBUG1 |
6 |
DEBUG |
5 |
INFO |
4 |
NOTE |
3 |
WARN |
2 |
ERROR |
1 |
QUIET |
0 |
Large value means more noise in the output file.
Note
Error and warning messages are written to stderr.
Each Logger object has its own output destination and verbose level. So
multiple Logger objects can be created to manage the message system without
affecting each other.
The methods provided by Logger class has the direct connection to the log level.
E.g. info()
print messages if the verbose level >= 4 (INFO):
>>> import sys
>>> from pyscf import lib
>>> log = lib.logger.Logger(sys.stdout, 4)
>>> log.info('info level')
info level
>>> log.verbose = 3
>>> log.info('info level')
>>> log.note('note level')
note level
timer¶
Logger object provides timer method for timing. Set TIMER_LEVEL
to
control at which level the timing information should be output. It is 5
(DEBUG) by default.
>>> import sys, time
>>> from pyscf import lib
>>> log = lib.logger.Logger(sys.stdout, 4)
>>> t0 = time.clock()
>>> log.timer('test', t0)
>>> lib.logger.TIMER_LEVEL = 4
>>> log.timer('test', t0)
CPU time for test 0.00 sec
-
class
pyscf.lib.logger.
Logger
(stdout=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, verbose=3)¶ Bases:
object
- Attributes:
- stdoutfile object or sys.stdout
The file to dump output message.
- verboseint
Large value means more noise in the output file.
-
debug
(msg, *args)¶
-
debug1
(msg, *args)¶
-
debug2
(msg, *args)¶
-
debug3
(msg, *args)¶
-
debug4
(msg, *args)¶
-
error
(msg, *args)¶
-
info
(msg, *args)¶
-
log
(msg, *args)¶
-
note
(msg, *args)¶
-
timer
(msg, cpu0=None, wall0=None)¶
-
timer_debug1
(msg, cpu0=None, wall0=None)¶
-
warn
(msg, *args)¶
-
pyscf.lib.logger.
debug
(rec, msg, *args)¶
-
pyscf.lib.logger.
debug1
(rec, msg, *args)¶
-
pyscf.lib.logger.
debug2
(rec, msg, *args)¶
-
pyscf.lib.logger.
debug3
(rec, msg, *args)¶
-
pyscf.lib.logger.
debug4
(rec, msg, *args)¶
-
pyscf.lib.logger.
error
(rec, msg, *args)¶
-
pyscf.lib.logger.
flush
(rec, msg, *args)¶
-
pyscf.lib.logger.
info
(rec, msg, *args)¶
-
pyscf.lib.logger.
log
(rec, msg, *args)¶
-
pyscf.lib.logger.
new_logger
(rec=None, verbose=None)¶ Create and return a
Logger
object- Args:
rec : An object which carries the attributes stdout and verbose
- verbosea Logger object, or integer or None
The verbose level. If verbose is a Logger object, the Logger object is returned. If verbose is not specified (None), rec.verbose will be used in the new Logger object.
-
pyscf.lib.logger.
note
(rec, msg, *args)¶
-
pyscf.lib.logger.
stdout
(rec, msg, *args)¶
-
pyscf.lib.logger.
timer
(rec, msg, cpu0=None, wall0=None)¶
-
pyscf.lib.logger.
timer_debug1
(rec, msg, cpu0=None, wall0=None)¶
-
pyscf.lib.logger.
warn
(rec, msg, *args)¶
pyscf.lib.misc module¶
Some helper functions
-
class
pyscf.lib.misc.
H5TmpFile
(filename=None, mode='a', *args, **kwargs)¶ Bases:
h5py._hl.files.File
Create and return an HDF5 temporary file.
- Kwargs:
- filenamestr or None
If a string is given, an HDF5 file of the given filename will be created. The temporary file will exist even if the H5TmpFile object is released. If nothing is specified, the HDF5 temporary file will be deleted when the H5TmpFile object is released.
The return object is an h5py.File object. The file will be automatically deleted when it is closed or the object is released (unless filename is specified).
Examples:
>>> from pyscf import lib >>> ftmp = lib.H5TmpFile()
-
exception
pyscf.lib.misc.
ProcessRuntimeError
¶ Bases:
RuntimeError
-
class
pyscf.lib.misc.
ProcessWithReturnValue
(group=None, target=None, name=None, args=(), kwargs=None)¶ Bases:
multiprocessing.context.Process
-
get
()¶ Wait until child process terminates
-
join
()¶ Wait until child process terminates
-
-
class
pyscf.lib.misc.
SinglePointScanner
¶ Bases:
object
-
class
pyscf.lib.misc.
StreamObject
¶ Bases:
object
For most methods, there are three stream functions to pipe computing stream:
1
.set_
function to update object attributes, egmf = scf.RHF(mol).set(conv_tol=1e-5)
is identical to proceed in two stepsmf = scf.RHF(mol); mf.conv_tol=1e-5
2
.run
function to execute the kenerl function (the function arguments are passed to kernel function). If keyword arguments is given, it will first call.set
function to update object attributes then execute the kernel function. Egmf = scf.RHF(mol).run(dm_init, conv_tol=1e-5)
is identical to three stepsmf = scf.RHF(mol); mf.conv_tol=1e-5; mf.kernel(dm_init)
3
.apply
function to apply the given function/class to the current object (function arguments and keyword arguments are passed to the given function). Egmol.apply(scf.RHF).run().apply(mcscf.CASSCF, 6, 4, frozen=4)
is identical tomf = scf.RHF(mol); mf.kernel(); mcscf.CASSCF(mf, 6, 4, frozen=4)
-
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.
-
check_sanity
()¶ Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.
-
kernel
(*args, **kwargs)¶ Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).
-
post_kernel
(envs)¶ A hook to be run after the main body of the kernel function. Internal variables are exposed to post_kernel through the “envs” dictionary. Return value of post_kernel function is not required.
-
pre_kernel
(envs)¶ A hook to be run before the main body of kernel function is executed. Internal variables are exposed to pre_kernel through the “envs” dictionary. Return value of pre_kernel function is not required.
-
run
(*args, **kwargs)¶ Call the kernel function of current object. args will be passed to kernel function. kwargs will be used to update the attributes of current object. The return value of method run is the object itself. This allows a series of functions/methods to be executed in pipe.
-
set
(*args, **kwargs)¶ Update the attributes of the current object. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.
-
stdout
= <_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>¶
-
verbose
= 0¶
-
view
(cls)¶ New view of object with the same attributes.
-
-
exception
pyscf.lib.misc.
ThreadRuntimeError
¶ Bases:
RuntimeError
-
class
pyscf.lib.misc.
ThreadWithReturnValue
(group=None, target=None, name=None, args=(), kwargs=None)¶ Bases:
threading.Thread
-
get
()¶ Wait until the thread terminates.
This blocks the calling thread until the thread whose join() method is called terminates – either normally or through an unhandled exception or until the optional timeout occurs.
When the timeout argument is present and not None, it should be a floating point number specifying a timeout for the operation in seconds (or fractions thereof). As join() always returns None, you must call is_alive() after join() to decide whether a timeout happened – if the thread is still alive, the join() call timed out.
When the timeout argument is not present or None, the operation will block until the thread terminates.
A thread can be join()ed many times.
join() raises a RuntimeError if an attempt is made to join the current thread as that would cause a deadlock. It is also an error to join() a thread before it has been started and attempts to do so raises the same exception.
-
join
()¶ Wait until the thread terminates.
This blocks the calling thread until the thread whose join() method is called terminates – either normally or through an unhandled exception or until the optional timeout occurs.
When the timeout argument is present and not None, it should be a floating point number specifying a timeout for the operation in seconds (or fractions thereof). As join() always returns None, you must call is_alive() after join() to decide whether a timeout happened – if the thread is still alive, the join() call timed out.
When the timeout argument is not present or None, the operation will block until the thread terminates.
A thread can be join()ed many times.
join() raises a RuntimeError if an attempt is made to join the current thread as that would cause a deadlock. It is also an error to join() a thread before it has been started and attempts to do so raises the same exception.
-
-
class
pyscf.lib.misc.
ThreadWithTraceBack
(group=None, target=None, name=None, args=(), kwargs=None)¶ Bases:
threading.Thread
-
join
()¶ Wait until the thread terminates.
This blocks the calling thread until the thread whose join() method is called terminates – either normally or through an unhandled exception or until the optional timeout occurs.
When the timeout argument is present and not None, it should be a floating point number specifying a timeout for the operation in seconds (or fractions thereof). As join() always returns None, you must call is_alive() after join() to decide whether a timeout happened – if the thread is still alive, the join() call timed out.
When the timeout argument is not present or None, the operation will block until the thread terminates.
A thread can be join()ed many times.
join() raises a RuntimeError if an attempt is made to join the current thread as that would cause a deadlock. It is also an error to join() a thread before it has been started and attempts to do so raises the same exception.
-
-
pyscf.lib.misc.
alias
(fn, alias_name=None)¶ The statement “fn1 = alias(fn)” in a class is equivalent to define the following method in the class:
Using alias function instead of fn1 = fn because some methods may be overloaded in the child class. Using “alias” can make sure that the overloaded mehods were called when calling the aliased method.
-
pyscf.lib.misc.
arg_first_match
(test, lst)¶
-
pyscf.lib.misc.
background
(func, *args, **kwargs)¶ applying function in background
-
pyscf.lib.misc.
background_process
(func, *args, **kwargs)¶ applying function in background
-
pyscf.lib.misc.
background_thread
(func, *args, **kwargs)¶ applying function in background
-
pyscf.lib.misc.
bg
(func, *args, **kwargs)¶ applying function in background
-
pyscf.lib.misc.
bg_process
(func, *args, **kwargs)¶ applying function in background
-
pyscf.lib.misc.
bg_thread
(func, *args, **kwargs)¶ applying function in background
-
pyscf.lib.misc.
bp
(func, *args, **kwargs)¶ applying function in background
-
pyscf.lib.misc.
c_double_arr
(m)¶
-
pyscf.lib.misc.
c_double_p
¶ alias of
pyscf.lib.misc.LP_c_double
-
pyscf.lib.misc.
c_int_arr
(m)¶
-
pyscf.lib.misc.
c_int_p
¶ alias of
pyscf.lib.misc.LP_c_int
-
pyscf.lib.misc.
c_null_ptr
¶ alias of
pyscf.lib.misc.LP_c_void_p
-
class
pyscf.lib.misc.
call_in_background
(*fns, **kwargs)¶ Bases:
object
Within this macro, function(s) can be executed asynchronously (the given functions are executed in background).
- Attributes:
- sync (bool): Whether to run in synchronized mode. The default value
is False (asynchoronized mode).
Examples:
>>> with call_in_background(fun) as async_fun: ... async_fun(a, b) # == fun(a, b) ... do_something_else()
>>> with call_in_background(fun1, fun2) as (afun1, afun2): ... afun2(a, b) ... do_something_else() ... afun2(a, b) ... do_something_else() ... afun1(a, b) ... do_something_else()
-
class
pyscf.lib.misc.
capture_stdout
¶ Bases:
object
redirect all stdout (c printf & python print) into a string
Examples:
>>> import os >>> from pyscf import lib >>> with lib.capture_stdout() as out: ... os.system('ls') >>> print(out.read())
-
read
()¶
-
-
pyscf.lib.misc.
check_sanity
(obj, keysref, stdout=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>)¶ Check misinput of class attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”.
-
pyscf.lib.misc.
class_as_method
(cls)¶ The statement “fn1 = alias(Class)” is equivalent to:
-
pyscf.lib.misc.
ctypes_stdout
¶ alias of
pyscf.lib.misc.capture_stdout
-
pyscf.lib.misc.
current_memory
()¶ Return the size of used memory and allocated virtual memory (in MB)
-
pyscf.lib.misc.
f_double_arr
(m)¶
-
pyscf.lib.misc.
f_int_arr
(m)¶
-
pyscf.lib.misc.
find_if
(test, lst)¶
-
pyscf.lib.misc.
finger
(a)¶ Fingerprint of numpy array
-
pyscf.lib.misc.
fingerprint
(a)¶ Fingerprint of numpy array
-
pyscf.lib.misc.
flatten
(lst)¶ flatten nested lists x[0] + x[1] + x[2] + …
Examples:
>>> flatten([[0, 2], [1], [[9, 8, 7]]]) [0, 2, 1, [9, 8, 7]]
-
pyscf.lib.misc.
fp
(a)¶ Fingerprint of numpy array
-
pyscf.lib.misc.
index_tril_to_pair
(ij)¶ Given tril-index ij, compute the pair indices (i,j) which satisfy ij = i * (i+1) / 2 + j
-
pyscf.lib.misc.
izip
(*args)¶ python2 izip == python3 zip
-
class
pyscf.lib.misc.
light_speed
(c)¶ Bases:
pyscf.lib.misc.temporary_env
Within the context of this macro, the environment varialbe LIGHT_SPEED can be customized.
Examples:
>>> with light_speed(15.): ... print(lib.param.LIGHT_SPEED) 15. >>> print(lib.param.LIGHT_SPEED) 137.03599967994
-
pyscf.lib.misc.
load_library
(libname)¶
-
pyscf.lib.misc.
member
(test, x, lst)¶
-
pyscf.lib.misc.
ndpointer
(*args, **kwargs)¶
-
pyscf.lib.misc.
num_threads
(n=None)¶ Set the number of OMP threads. If argument is not specified, the function will return the total number of available OMP threads.
It’s recommended to call this function to set OMP threads than “os.environ[‘OMP_NUM_THREADS’] = int(n)”. This is because environment variables like OMP_NUM_THREADS were read when a module was imported. They cannot be reset through os.environ after the module was loaded.
Examples:
>>> from pyscf import lib >>> print(lib.num_threads()) 8 >>> lib.num_threads(4) 4 >>> print(lib.num_threads()) 4
-
class
pyscf.lib.misc.
omnimethod
(func)¶ Bases:
object
-
pyscf.lib.misc.
overwrite_mro
(obj, mro)¶ A hacky function to overwrite the __mro__ attribute
-
pyscf.lib.misc.
prange
(start, end, step)¶ This function splits the number sequence between “start” and “end” using uniform “step” length. It yields the boundary (start, end) for each fragment.
Examples:
>>> for p0, p1 in lib.prange(0, 8, 2): ... print(p0, p1) (0, 2) (2, 4) (4, 6) (6, 8)
-
pyscf.lib.misc.
prange_tril
(start, stop, blocksize)¶ Similar to
prange()
, yeilds start (p0) and end (p1) with the restriction p1*(p1+1)/2-p0*(p0+1)/2 < blocksizeExamples:
>>> for p0, p1 in lib.prange_tril(0, 10, 25): ... print(p0, p1) (0, 6) (6, 9) (9, 10)
-
class
pyscf.lib.misc.
quite_run
¶ Bases:
object
capture all stdout (c printf & python print) but output nothing
Examples:
>>> import os >>> from pyscf import lib >>> with lib.quite_run(): ... os.system('ls')
-
pyscf.lib.misc.
remove_dup
(test, lst, from_end=False)¶
-
pyscf.lib.misc.
remove_if
(test, lst)¶
-
pyscf.lib.misc.
square_mat_in_trilu_indices
(n)¶ Return a n x n symmetric index matrix, in which the elements are the indices of the unique elements of a tril vector [0 1 3 … ] [1 2 4 … ] [3 4 5 … ] [… ]
-
class
pyscf.lib.misc.
temporary_env
(obj, **kwargs)¶ Bases:
object
Within the context of this macro, the attributes of the object are temporarily updated. When the program goes out of the scope of the context, the original value of each attribute will be restored.
Examples:
>>> with temporary_env(lib.param, LIGHT_SPEED=15., BOHR=2.5): ... print(lib.param.LIGHT_SPEED, lib.param.BOHR) 15. 2.5 >>> print(lib.param.LIGHT_SPEED, lib.param.BOHR) 137.03599967994 0.52917721092
-
pyscf.lib.misc.
tril_product
(*iterables, **kwds)¶ Cartesian product in lower-triangular form for multiple indices
For a given list of indices (iterables), this function yields all indices such that the sub-indices given by the kwarg tril_idx satisfy a lower-triangular form. The lower-triangular form satisfies:
\[i[tril_idx[0]] >= i[tril_idx[1]] >= ... >= i[tril_idx[len(tril_idx)-1]]\]- Args:
*iterables: Variable length argument list of indices for the cartesian product **kwds: Arbitrary keyword arguments. Acceptable keywords include:
repeat (int): Number of times to repeat the iterables tril_idx (array_like): Indices to put into lower-triangular form.
- Yields:
product (tuple): Tuple in lower-triangular form.
- Examples:
Specifying no tril_idx is equivalent to just a cartesian product.
>>> list(tril_product(range(2), repeat=2)) [(0, 0), (0, 1), (1, 0), (1, 1)]
We can specify only sub-indices to satisfy a lower-triangular form:
>>> list(tril_product(range(2), repeat=3, tril_idx=[1,2])) [(0, 0, 0), (0, 1, 0), (0, 1, 1), (1, 0, 0), (1, 1, 0), (1, 1, 1)]
We specify all indices to satisfy a lower-triangular form, useful for iterating over the symmetry unique elements of occupied/virtual orbitals in a 3-particle operator:
>>> list(tril_product(range(3), repeat=3, tril_idx=[0,1,2])) [(0, 0, 0), (1, 0, 0), (1, 1, 0), (1, 1, 1), (2, 0, 0), (2, 1, 0), (2, 1, 1), (2, 2, 0), (2, 2, 1), (2, 2, 2)]
-
pyscf.lib.misc.
with_doc
(doc)¶ Use this decorator to add doc string for function
@with_doc(doc) def fn:
…
is equivalent to
fn.__doc__ = doc
-
class
pyscf.lib.misc.
with_omp_threads
(nthreads=None)¶ Bases:
object
Using this macro to create a temporary context in which the number of OpenMP threads are set to the required value. When the program exits the context, the number OpenMP threads will be restored.
- Args:
nthreads : int
Examples:
>>> from pyscf import lib >>> print(lib.num_threads()) 8 >>> with lib.with_omp_threads(2): ... print(lib.num_threads()) 2 >>> print(lib.num_threads()) 8
pyscf.lib.numpy_helper module¶
Extension to numpy and scipy
-
class
pyscf.lib.numpy_helper.
NPArrayWithTag
¶ Bases:
numpy.ndarray
-
pyscf.lib.numpy_helper.
cartesian_prod
(arrays, out=None)¶ Generate a cartesian product of input arrays. http://stackoverflow.com/questions/1208118/using-numpy-to-build-an-array-of-all-combinations-of-two-arrays
- Args:
- arrayslist of array-like
1-D arrays to form the cartesian product of.
- outndarray
Array to place the cartesian product in.
- Returns:
- outndarray
2-D array of shape (M, len(arrays)) containing cartesian products formed of input arrays.
Examples:
>>> cartesian_prod(([1, 2, 3], [4, 5], [6, 7])) array([[1, 4, 6], [1, 4, 7], [1, 5, 6], [1, 5, 7], [2, 4, 6], [2, 4, 7], [2, 5, 6], [2, 5, 7], [3, 4, 6], [3, 4, 7], [3, 5, 6], [3, 5, 7]])
-
pyscf.lib.numpy_helper.
cond
(x, p=None)¶ Compute the condition number
-
pyscf.lib.numpy_helper.
condense
(opname, a, locs)¶ nd = loc[-1] out = numpy.empty((nd,nd)) for i,i0 in enumerate(loc): i1 = loc[i+1] for j,j0 in enumerate(loc): j1 = loc[j+1] out[i,j] = op(a[i0:i1,j0:j1]) return out
-
pyscf.lib.numpy_helper.
ddot
(a, b, alpha=1, c=None, beta=0)¶ Matrix-matrix multiplication for double precision arrays
-
pyscf.lib.numpy_helper.
direct_sum
(subscripts, *operands)¶ Apply the summation over many operands with the einsum fashion.
Examples:
>>> a = numpy.random.random((6,5)) >>> b = numpy.random.random((4,3,2)) >>> direct_sum('ij,klm->ijklm', a, b).shape (6, 5, 4, 3, 2) >>> direct_sum('ij,klm', a, b).shape (6, 5, 4, 3, 2) >>> direct_sum('i,j,klm->mjlik', a[0], a[:,0], b).shape (2, 6, 3, 5, 4) >>> direct_sum('ij-klm->ijklm', a, b).shape (6, 5, 4, 3, 2) >>> direct_sum('ij+klm', a, b).shape (6, 5, 4, 3, 2) >>> direct_sum('-i-j+klm->mjlik', a[0], a[:,0], b).shape (2, 6, 3, 5, 4) >>> c = numpy.random((3,5)) >>> z = direct_sum('ik+jk->kij', a, c).shape # This is slow >>> abs(a.T.reshape(5,6,1) + c.reshape(5,1,3) - z).sum() 0.0
-
pyscf.lib.numpy_helper.
dot
(a, b, alpha=1, c=None, beta=0)¶
-
pyscf.lib.numpy_helper.
einsum
(subscripts, *tensors, **kwargs)¶ Perform a more efficient einsum via reshaping to a matrix multiply.
Current differences compared to numpy.einsum: This assumes that each repeated index is actually summed (i.e. no ‘i,i->i’) and appears only twice (i.e. no ‘ij,ik,il->jkl’). The output indices must be explicitly specified (i.e. ‘ij,j->i’ and not ‘ij,j’).
-
pyscf.lib.numpy_helper.
expm
(a)¶ Equivalent to scipy.linalg.expm
-
pyscf.lib.numpy_helper.
frompointer
(pointer, count, dtype=<class 'float'>)¶ Interpret a buffer that the pointer refers to as a 1-dimensional array.
- Args:
- pointerint or ctypes pointer
address of a buffer
- countint
Number of items to read.
- dtypedata-type, optional
Data-type of the returned array; default: float.
Examples:
>>> s = numpy.ones(3, dtype=numpy.int32) >>> ptr = s.ctypes.data >>> frompointer(ptr, count=6, dtype=numpy.int16) [1, 0, 1, 0, 1, 0]
-
pyscf.lib.numpy_helper.
hermi_sum
(a, axes=None, hermi=1, inplace=False, out=None)¶ Computing a + a.T.conj() with better memory efficiency
Examples:
>>> transpose_sum(numpy.arange(4.).reshape(2,2)) [[ 0. 3.] [ 3. 6.]]
-
pyscf.lib.numpy_helper.
hermi_triu
(mat, hermi=1, inplace=True)¶ Use the elements of the lower triangular part to fill the upper triangular part.
- Kwargs:
filltriu : int
1 (default) return a hermitian matrix2 return an anti-hermitian matrix
Examples:
>>> unpack_row(numpy.arange(9.).reshape(3,3), 1) [[ 0. 3. 6.] [ 3. 4. 7.] [ 6. 7. 8.]] >>> unpack_row(numpy.arange(9.).reshape(3,3), 2) [[ 0. -3. -6.] [ 3. 4. -7.] [ 6. 7. 8.]]
-
pyscf.lib.numpy_helper.
pack_tril
(mat, axis=- 1, out=None)¶ flatten the lower triangular part of a matrix. Given mat, it returns mat[…,numpy.tril_indices(mat.shape[0])]
Examples:
>>> pack_tril(numpy.arange(9).reshape(3,3)) [0 3 4 6 7 8]
-
pyscf.lib.numpy_helper.
solve_lineq_by_SVD
(a, b)¶ Solving a * x = b. If a is a singular matrix, its small SVD values are neglected.
-
pyscf.lib.numpy_helper.
split_reshape
(a, shapes)¶ Split a vector into multiple tensors. shapes is a list of tuples. The entries of shapes indicate the shape of each tensor.
- Returns:
tensors : a list of tensors
Examples:
>>> a = numpy.arange(12) >>> split_reshape(a, ((2,3), (1,), ((2,2), (1,1)))) [array([[0, 1, 2], [3, 4, 5]]), array([6]), [array([[ 7, 8], [ 9, 10]]), array([[11]])]]
-
pyscf.lib.numpy_helper.
tag_array
(a, **kwargs)¶ Attach attributes to numpy ndarray. The attribute name and value are obtained from the keyword arguments.
-
pyscf.lib.numpy_helper.
take_2d
(a, idx, idy, out=None)¶ Equivalent to a[idx[:,None],idy] for a 2D array.
Examples:
>>> out = numpy.arange(9.).reshape(3,3) >>> take_2d(a, [0,2], [0,2]) [[ 0. 2.] [ 6. 8.]]
-
pyscf.lib.numpy_helper.
takebak_2d
(out, a, idx, idy)¶ Reverse operation of take_2d. Equivalent to out[idx[:,None],idy] += a for a 2D array.
Examples:
>>> out = numpy.zeros((3,3)) >>> takebak_2d(out, numpy.ones((2,2)), [0,2], [0,2]) [[ 1. 0. 1.] [ 0. 0. 0.] [ 1. 0. 1.]]
-
pyscf.lib.numpy_helper.
transpose
(a, axes=None, inplace=False, out=None)¶ Transposing an array with better memory efficiency
Examples:
>>> transpose(numpy.ones((3,2))) [[ 1. 1. 1.] [ 1. 1. 1.]]
-
pyscf.lib.numpy_helper.
transpose_sum
(a, inplace=False, out=None)¶ Computing a + a.T with better memory efficiency
Examples:
>>> transpose_sum(numpy.arange(4.).reshape(2,2)) [[ 0. 3.] [ 3. 6.]]
-
pyscf.lib.numpy_helper.
unpack_row
(tril, row_id)¶ Extract one row of the lower triangular part of a matrix. It is equivalent to unpack_tril(a)[row_id]
Examples:
>>> unpack_row(numpy.arange(6.), 0) [ 0. 1. 3.] >>> unpack_tril(numpy.arange(6.))[0] [ 0. 1. 3.]
-
pyscf.lib.numpy_helper.
unpack_tril
(tril, filltriu=1, axis=- 1, out=None)¶ Reversed operation of pack_tril.
- Kwargs:
filltriu : int
0 Do not fill the upper triangular part, random number may appear in the upper triangular part1 (default) Transpose the lower triangular part to fill the upper triangular part2 Similar to filltriu=1, negative of the lower triangular part is assign to the upper triangular part to make the matrix anti-hermitian
Examples:
>>> unpack_tril(numpy.arange(6.)) [[ 0. 1. 3.] [ 1. 2. 4.] [ 3. 4. 5.]] >>> unpack_tril(numpy.arange(6.), 0) [[ 0. 0. 0.] [ 1. 2. 0.] [ 3. 4. 5.]] >>> unpack_tril(numpy.arange(6.), 2) [[ 0. -1. -3.] [ 1. 2. -4.] [ 3. 4. 5.]]
-
pyscf.lib.numpy_helper.
zdot
(a, b, alpha=1, c=None, beta=0)¶ Matrix-matrix multiplication for double complex arrays
pyscf.lib.parameters module¶
PySCF environment variables are defined in this module.
Scratch directory¶
The PySCF scratch directory is specified by TMPDIR
. Its default value
is the same to the system-wide environment variable TMPDIR
. It can be
overwritten by the environment variable PYSCF_TMPDIR
. Another place to set
TMPDIR
is the global configuration file (see global_config).
Maximum memory¶
The variable MAX_MEMORY
defines the maximum memory that PySCF can be
used in the calculation. Its unit is MB. The default value is 4000 MB. It can
be overwritten by the system-wide environment variable PYSCF_MAX_MEMORY
.
MAX_MEMORY
can also be set in global_config file.
Note
Some calculations may exceed the max_memory limit, especially
when the attribute Mole.incore_anyway
was set.
pyscf.lib.tblis_einsum module¶
Module contents¶
C extensions and helper functions