pyscf.tools package#
Submodules#
pyscf.tools.c60struct module#
pyscf.tools.chgcar module#
Vasp CHGCAR file format
See also https://cms.mpi.univie.ac.at/vasp/vasp/CHGCAR_file.html
- class pyscf.tools.chgcar.CHGCAR(cell, nx=60, ny=60, nz=60, resolution=None, margin=3.0)[source]#
Bases:
Cube
Read-write of the Vasp CHGCAR files
- pyscf.tools.chgcar.density(cell, outfile, dm, nx=60, ny=60, nz=60, resolution=None)[source]#
Calculates electron density and write out in CHGCAR format.
- Args:
- cellMole or Cell object
Mole or pbc Cell. If Mole object is given, the program will guess a cubic lattice for the molecule.
- outfilestr
Name of Cube file to be written.
- dmndarray
Density matrix of molecule.
- Kwargs:
- nxint
Number of grid point divisions in x direction. Note this is function of the molecule’s size; a larger molecule will have a coarser representation than a smaller one for the same value.
- nyint
Number of grid point divisions in y direction.
- nzint
Number of grid point divisions in z direction.
- Returns:
No return value. This function outputs a VASP chgcarlike file (with phase if desired)…it can be opened in VESTA or VMD or many other softwares
Examples:
>>> # generates the first MO from the list of mo_coefficents >>> from pyscf.pbc import gto, scf >>> from pyscf.tools import chgcar >>> cell = gto.M(atom='H 0 0 0; H 0 0 1', a=numpy.eye(3)*3) >>> mf = scf.RHF(cell).run() >>> chgcar.density(cell, 'h2.CHGCAR', mf.make_rdm1())
- pyscf.tools.chgcar.orbital(cell, outfile, coeff, nx=60, ny=60, nz=60, resolution=None)[source]#
Calculate orbital value on real space grid and write out in CHGCAR format.
- Args:
- cellMole or Cell object
Mole or pbc Cell. If Mole object is given, the program will guess a cubic lattice for the molecule.
- outfilestr
Name of Cube file to be written.
- dmndarray
Density matrix of molecule.
- Kwargs:
- nxint
Number of grid point divisions in x direction. Note this is function of the molecule’s size; a larger molecule will have a coarser representation than a smaller one for the same value.
- nyint
Number of grid point divisions in y direction.
- nzint
Number of grid point divisions in z direction.
- Returns:
No return value. This function outputs a VASP chgcarlike file (with phase if desired)…it can be opened in VESTA or VMD or many other softwares
Examples:
>>> # generates the first MO from the list of mo_coefficents >>> from pyscf.pbc import gto, scf >>> from pyscf.tools import chgcar >>> cell = gto.M(atom='H 0 0 0; H 0 0 1', a=numpy.eye(3)*3) >>> mf = scf.RHF(cell).run() >>> chgcar.orbital(cell, 'h2_mo1.CHGCAR', mf.mo_coeff[:,0])
pyscf.tools.chkfile_util module#
- pyscf.tools.chkfile_util.dump_mo(filename, key='scf')[source]#
Read scf/mcscf information from chkfile, then dump the orbital coefficients.
pyscf.tools.cubegen module#
Gaussian cube file format. Reference: http://paulbourke.net/dataformats/cube/ https://h5cube-spec.readthedocs.io/en/latest/cubeformat.html http://gaussian.com/cubegen/
The output cube file has the following format
Comment line Comment line N_atom Ox Oy Oz # number of atoms, followed by the coordinates of the origin N1 vx1 vy1 vz1 # number of grids along each axis, followed by the step size in x/y/z direction. N2 vx2 vy2 vz2 # … N3 vx3 vy3 vz3 # … Atom1 Z1 x y z # Atomic number, charge, and coordinates of the atom … # … AtomN ZN x y z # … Data on grids # (N1*N2) lines of records, each line has N3 elements
- class pyscf.tools.cubegen.Cube(mol, nx=80, ny=80, nz=80, resolution=None, margin=3.0, origin=None, extent=None)[source]#
Bases:
object
Read-write of the Gaussian CUBE files
- Attributes:
- nxint
Number of grid point divisions in x direction. Note this is function of the molecule’s size; a larger molecule will have a coarser representation than a smaller one for the same value. Conflicts to keyword resolution.
- nyint
Number of grid point divisions in y direction.
- nzint
Number of grid point divisions in z direction.
- resolution: float
Resolution of the mesh grid in the cube box. If resolution is given in the input, the input nx/ny/nz have no effects. The value of nx/ny/nz will be determined by the resolution and the cube box size. The unit is Bohr.
- pyscf.tools.cubegen.density(mol, outfile, dm, nx=80, ny=80, nz=80, resolution=None, margin=3.0)[source]#
Calculates electron density and write out in cube format.
- Args:
- molMole
Molecule to calculate the electron density for.
- outfilestr
Name of Cube file to be written.
- dmndarray
Density matrix of molecule.
- Kwargs:
- nxint
Number of grid point divisions in x direction. Note this is function of the molecule’s size; a larger molecule will have a coarser representation than a smaller one for the same value. Conflicts to keyword resolution.
- nyint
Number of grid point divisions in y direction.
- nzint
Number of grid point divisions in z direction.
- resolution: float
Resolution of the mesh grid in the cube box. If resolution is given in the input, the input nx/ny/nz have no effects. The value of nx/ny/nz will be determined by the resolution and the cube box size.
- pyscf.tools.cubegen.mep(mol, outfile, dm, nx=80, ny=80, nz=80, resolution=None, margin=3.0)[source]#
Calculates the molecular electrostatic potential (MEP) and write out in cube format.
- Args:
- molMole
Molecule to calculate the electron density for.
- outfilestr
Name of Cube file to be written.
- dmndarray
Density matrix of molecule.
- Kwargs:
- nxint
Number of grid point divisions in x direction. Note this is function of the molecule’s size; a larger molecule will have a coarser representation than a smaller one for the same value. Conflicts to keyword resolution.
- nyint
Number of grid point divisions in y direction.
- nzint
Number of grid point divisions in z direction.
- resolution: float
Resolution of the mesh grid in the cube box. If resolution is given in the input, the input nx/ny/nz have no effects. The value of nx/ny/nz will be determined by the resolution and the cube box size.
- pyscf.tools.cubegen.orbital(mol, outfile, coeff, nx=80, ny=80, nz=80, resolution=None, margin=3.0)[source]#
Calculate orbital value on real space grid and write out in cube format.
- Args:
- molMole
Molecule to calculate the electron density for.
- outfilestr
Name of Cube file to be written.
- coeff1D array
coeff coefficient.
- Kwargs:
- nxint
Number of grid point divisions in x direction. Note this is function of the molecule’s size; a larger molecule will have a coarser representation than a smaller one for the same value. Conflicts to keyword resolution.
- nyint
Number of grid point divisions in y direction.
- nzint
Number of grid point divisions in z direction.
- resolution: float
Resolution of the mesh grid in the cube box. If resolution is given in the input, the input nx/ny/nz have no effects. The value of nx/ny/nz will be determined by the resolution and the cube box size.
pyscf.tools.dump_mat module#
- pyscf.tools.dump_mat.dump_mo(mol, c, label=None, ncol=5, digits=5, start=0)[source]#
Format print for orbitals
- Args:
- stdoutfile object
eg sys.stdout, or stdout = open(‘/path/to/file’) or mol.stdout if mol is an object initialized from
gto.Mole
- cnumpy.ndarray
Orbitals, each column is an orbital
- Kwargs:
- labellist of strings
Row labels (default is AO labels)
Examples:
>>> from pyscf import gto >>> mol = gto.M(atom='C 0 0 0') >>> mo = numpy.eye(mol.nao_nr()) >>> dump_mo(mol, mo) #0 #1 #2 #3 #4 #5 #6 #7 #8 0 C 1s 1.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0 C 2s 0.00 1.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0 C 3s 0.00 0.00 1.00 0.00 0.00 0.00 0.00 0.00 0.00 0 C 2px 0.00 0.00 0.00 1.00 0.00 0.00 0.00 0.00 0.00 0 C 2py 0.00 0.00 0.00 0.00 1.00 0.00 0.00 0.00 0.00 0 C 2pz 0.00 0.00 0.00 0.00 0.00 1.00 0.00 0.00 0.00 0 C 3px 0.00 0.00 0.00 0.00 0.00 0.00 1.00 0.00 0.00 0 C 3py 0.00 0.00 0.00 0.00 0.00 0.00 0.00 1.00 0.00 0 C 3pz 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 1.00
- pyscf.tools.dump_mat.dump_rec(stdout, c, label=None, label2=None, ncol=5, digits=5, start=0)[source]#
Print an array in rectangular format
- Args:
- stdoutfile object
eg sys.stdout, or stdout = open(‘/path/to/file’) or mol.stdout if mol is an object initialized from
gto.Mole
- cnumpy.ndarray
coefficients
- Kwargs:
- labellist of strings
Row labels (default is 1,2,3,4,…)
- label2list of strings
Col labels (default is 1,2,3,4,…)
- ncolint
Number of columns in the format output (default 5)
- digitsint
Number of digits of precision for floating point output (default 5)
- startint
The number to start to count the index (default 0)
Examples:
>>> import sys, numpy >>> dm = numpy.eye(3) >>> dump_rec(sys.stdout, dm) #0 #1 #2 0 1.00000 0.00000 0.00000 1 0.00000 1.00000 0.00000 2 0.00000 0.00000 1.00000 >>> from pyscf import gto >>> mol = gto.M(atom='C 0 0 0') >>> dm = numpy.eye(mol.nao_nr()) >>> dump_rec(sys.stdout, dm, label=mol.ao_labels(), ncol=9, digits=2) #0 #1 #2 #3 #4 #5 #6 #7 #8 0 C 1s 1.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0 C 2s 0.00 1.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0 C 3s 0.00 0.00 1.00 0.00 0.00 0.00 0.00 0.00 0.00 0 C 2px 0.00 0.00 0.00 1.00 0.00 0.00 0.00 0.00 0.00 0 C 2py 0.00 0.00 0.00 0.00 1.00 0.00 0.00 0.00 0.00 0 C 2pz 0.00 0.00 0.00 0.00 0.00 1.00 0.00 0.00 0.00 0 C 3px 0.00 0.00 0.00 0.00 0.00 0.00 1.00 0.00 0.00 0 C 3py 0.00 0.00 0.00 0.00 0.00 0.00 0.00 1.00 0.00 0 C 3pz 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 1.00
- pyscf.tools.dump_mat.dump_tri(stdout, c, label=None, ncol=5, digits=5, start=0)[source]#
Format print for the lower triangular part of an array
- Args:
- stdoutfile object
eg sys.stdout, or stdout = open(‘/path/to/file’) or mol.stdout if mol is an object initialized from
gto.Mole
- cnumpy.ndarray
coefficients
- Kwargs:
- labellist of strings
Row labels (default is 1,2,3,4,…)
- ncolint
Number of columns in the format output (default 5)
- digitsint
Number of digits of precision for floating point output (default 5)
- startint
The number to start to count the index (default 0)
Examples:
>>> import sys, numpy >>> dm = numpy.eye(3) >>> dump_tri(sys.stdout, dm) #0 #1 #2 0 1.00000 1 0.00000 1.00000 2 0.00000 0.00000 1.00000 >>> from pyscf import gto >>> mol = gto.M(atom='C 0 0 0') >>> dm = numpy.eye(mol.nao_nr()) >>> dump_tri(sys.stdout, dm, label=mol.ao_labels(), ncol=9, digits=2) #0 #1 #2 #3 #4 #5 #6 #7 #8 0 C 1s 1.00 0 C 2s 0.00 1.00 0 C 3s 0.00 0.00 1.00 0 C 2px 0.00 0.00 0.00 1.00 0 C 2py 0.00 0.00 0.00 0.00 1.00 0 C 2pz 0.00 0.00 0.00 0.00 0.00 1.00 0 C 3px 0.00 0.00 0.00 0.00 0.00 0.00 1.00 0 C 3py 0.00 0.00 0.00 0.00 0.00 0.00 0.00 1.00 0 C 3pz 0.00 0.00 0.00 0.00 0.00 0.00 0.00 0.00 1.00
pyscf.tools.fcidump module#
FCIDUMP functions (write, read) for real Hamiltonian
- pyscf.tools.fcidump.from_chkfile(filename, chkfile, tol=1e-15, float_format=' %.16g', molpro_orbsym=False, orbsym=None)[source]#
Read SCF results from PySCF chkfile and transform 1-electron, 2-electron integrals using the SCF orbitals. The transformed integrals is written to FCIDUMP
- Kwargs:
- molpro_orbsym (bool): Whether to dump the orbsym in Molpro orbsym
convention as documented in https://www.molpro.net/info/current/doc/manual/node36.html
- pyscf.tools.fcidump.from_integrals(filename, h1e, h2e, nmo, nelec, nuc=0, ms=0, orbsym=None, tol=1e-15, float_format=' %.16g')[source]#
Convert the given 1-electron and 2-electron integrals to FCIDUMP format
- pyscf.tools.fcidump.from_mcscf(mc, filename, tol=1e-15, float_format=' %.16g', molpro_orbsym=False)[source]#
Use the given MCSCF object to obtain the CAS 1-electron and 2-electron integrals and dump them to FCIDUMP.
- Kwargs:
tol (float): Threshold for writing elements to FCIDUMP float_format (str): Float format for writing elements to FCIDUMP molpro_orbsym (bool): Whether to dump the orbsym in Molpro orbsym
convention as documented in https://www.molpro.net/manual/doku.php?id=general_program_structure#symmetry
- pyscf.tools.fcidump.from_mo(mol, filename, mo_coeff, orbsym=None, tol=1e-15, float_format=' %.16g', molpro_orbsym=False, ms=0)[source]#
Use the given MOs to transform the 1-electron and 2-electron integrals then dump them to FCIDUMP.
- Kwargs:
- molpro_orbsym (bool): Whether to dump the orbsym in Molpro orbsym
convention as documented in https://www.molpro.net/info/current/doc/manual/node36.html
- pyscf.tools.fcidump.from_scf(mf, filename, tol=1e-15, float_format=' %.16g', molpro_orbsym=False)[source]#
Use the given SCF object to transform the 1-electron and 2-electron integrals then dump them to FCIDUMP.
- Kwargs:
- molpro_orbsym (bool): Whether to dump the orbsym in Molpro orbsym
convention as documented in https://www.molpro.net/info/current/doc/manual/node36.html
- pyscf.tools.fcidump.read(filename, molpro_orbsym=False, verbose=True)[source]#
Parse FCIDUMP. Return a dictionary to hold the integrals and parameters with keys: H1, H2, ECORE, NORB, NELEC, MS, ORBSYM, ISYM
- Kwargs:
- molpro_orbsym (bool): Whether the orbsym in the FCIDUMP file is in
Molpro orbsym convention as documented in:
https://www.molpro.net/info/current/doc/manual/node36.html
In return, orbsym is converted to pyscf symmetry convention
verbose (bool): Whether to print debugging information
- pyscf.tools.fcidump.scf_from_fcidump(mf, filename, molpro_orbsym=False)[source]#
Update the SCF object with the quantities defined in FCIDUMP file
pyscf.tools.mo_mapping module#
- pyscf.tools.mo_mapping.mo_1to1map(s)[source]#
Given <i|j>, search for the 1-to-1 mapping between i and j.
- Returns:
a list [j-close-to-i for i in <bra|]
- pyscf.tools.mo_mapping.mo_comps(aolabels_or_baslst, mol, mo_coeff, cart=False, orth_method='meta_lowdin')[source]#
Given AO(s), show how the AO(s) are distributed in MOs.
- Args:
- aolabels_or_baslstfilter function or AO labels or AO index
If it’s a function, the AO indices are the items for which the function return value is true.
- Kwargs:
- cartbool
whether the orbital coefficients are based on cartesian basis.
- orth_methodstr
The localization method to generated orthogonal AO upon which the AO contribution are computed. It can be one of ‘meta_lowdin’, ‘lowdin’ or ‘nao’.
- Returns:
A list of float to indicate the total contributions (normalized to 1) of localized AOs
Examples:
>>> from pyscf import gto, scf >>> from pyscf.tools import mo_mapping >>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='6-31g') >>> mf = scf.RHF(mol).run() >>> comp = mo_mapping.mo_comps('F 2s', mol, mf.mo_coeff) >>> print('MO-id F-2s components') >>> for i,c in enumerate(comp): ... print('%-3d %.10f' % (i, c)) MO-id components 0 0.0000066344 1 0.8796915532 2 0.0590259826 3 0.0000000000 4 0.0000000000 5 0.0435028851 6 0.0155889103 7 0.0000000000 8 0.0000000000 9 0.0000822361 10 0.0021017982
pyscf.tools.molden module#
- pyscf.tools.molden.from_mo(mol, filename, mo_coeff, spin='Alpha', symm=None, ene=None, occ=None, ignore_h=True)[source]#
Dump the given MOs in Molden format
- pyscf.tools.molden.from_scf(mf, filename, ignore_h=True)[source]#
Dump the given SCF object in Molden format
- pyscf.tools.molden.orbital_coeff(mol, fout, mo_coeff, spin='Alpha', symm=None, ene=None, occ=None, ignore_h=True)[source]#
- pyscf.tools.molden.parse(moldenfile, verbose=0)#
Extract mol and orbitals from molden file
- pyscf.tools.molden.read(moldenfile, verbose=0)#
Extract mol and orbitals from molden file
- pyscf.tools.molden.remove_high_l(mol, mo_coeff=None)[source]#
Remove high angular momentum (l >= 5) functions before dumping molden file. If molden function raised error message
RuntimeError l=5 is not supported
, you can use this function to format orbitals.Note the formated orbitals may have normalization problem. Some visualization tool will complain about the orbital normalization error.
Examples:
>>> mol1, orb1 = remove_high_l(mol, mf.mo_coeff) >>> molden.from_mo(mol1, outputfile, orb1)
pyscf.tools.ring module#
pyscf.tools.wfn_format module#
GAMESS WFN File format