pyscf.symm package#
Submodules#
pyscf.symm.Dmatrix module#
Wigner rotation D-matrix for real spherical harmonics
- pyscf.symm.Dmatrix.Dmatrix(l, alpha, beta, gamma, reorder_p=False)[source]#
Wigner rotation D-matrix
D_{mm’} = <lm|R(alpha,beta,gamma)|lm’> alpha, beta, gamma are Euler angles (in z-y-z convention)
- Kwargs:
reorder_p (bool): Whether to put the p functions in the (x,y,z) order.
- pyscf.symm.Dmatrix.dmatrix(l, beta, reorder_p=False)[source]#
Wigner small-d matrix (in z-y-z convention)
- pyscf.symm.Dmatrix.get_euler_angles(c1, c2)[source]#
Find the three Euler angles (alpha, beta, gamma in z-y-z convention) that rotates coordinates c1 to coordinates c2.
yp = numpy.einsum(‘j,kj->k’, c1[1], geom.rotation_mat(c1[2], beta)) tmp = numpy.einsum(‘ij,kj->ik’, c1 , geom.rotation_mat(c1[2], alpha)) tmp = numpy.einsum(‘ij,kj->ik’, tmp, geom.rotation_mat(yp , beta )) c2 = numpy.einsum(‘ij,kj->ik’, tmp, geom.rotation_mat(c2[2], gamma))
(For backward compatibility) if c1 and c2 are two points in the real space, the Euler angles define the rotation transforms the old coordinates to the new coordinates (new_x, new_y, new_z) in which c1 is identical to c2.
tmp = numpy.einsum(‘j,kj->k’, c1 , geom.rotation_mat((0,0,1), gamma)) tmp = numpy.einsum(‘j,kj->k’, tmp, geom.rotation_mat((0,1,0), beta) ) c2 = numpy.einsum(‘j,kj->k’, tmp, geom.rotation_mat((0,0,1), alpha))
pyscf.symm.addons module#
- pyscf.symm.addons.eigh(h, orbsym)[source]#
Solve eigenvalue problem based on the symmetry information for basis. See also pyscf/lib/linalg_helper.py
eigh_by_blocks()
Examples:
>>> from pyscf import gto, 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) >>> symm.eigh(vnuc_so, orbsym) (array([-4.50766885, -1.80666351, -1.7808565 , -1.7808565 , -1.74189134, -0.98998583, -0.98998583, -0.40322226, -0.30242374, -0.07608981]), ...)
- pyscf.symm.addons.find_symmetric_mo(moso, ovlpso, thr=1e-08)[source]#
Find the list of MO of that thransform like particular irrep
- Args:
- moso2D float array
Overlap matrix of symmetry-adapted AO and MO, it can be obtained by reduce(numpy.dot, (csym.T.conj(), s, mo)), where csym taken from mol.symm_orb, and s is the AO overlap matrix
- ovlpso2D float array
Overlap matrix between symmetry-adapted AO, it can be obtained by reduce(numpy.dot, (csym.T.conj(), s, csym))
- Kwargs:
- thrfloat
Threshold to consider MO symmetry-adapted
- Returns:
1D bool array to select symmetry adapted MO
- pyscf.symm.addons.irrep_id2name(gpname, irrep_id)[source]#
Convert the internal irrep ID to irrep symbol
- Args:
- gpnamestr
The point group symbol
- irrep_idint
See IRREP_ID_TABLE in pyscf/symm/param.py
- Returns:
Irrep symbol, str
- pyscf.symm.addons.irrep_name2id(gpname, symb)[source]#
Convert the irrep symbol to internal irrep ID
- Args:
- gpnamestr
The point group symbol
- symbstr
Irrep symbol
- Returns:
Irrep ID, int
- pyscf.symm.addons.label_orb_symm(mol, irrep_name, symm_orb, mo, s=None, check=True, tol=1e-09)[source]#
Label the symmetry of given orbitals
irrep_name can be either the symbol or the ID of the irreducible representation. If the ID is provided, it returns the numeric code associated with XOR operator, see
symm.param.IRREP_ID_TABLE()
- Args:
mol : an instance of
Mole
- irrep_namelist of str or int
A list of irrep ID or name, it can be either mol.irrep_id or mol.irrep_name. It can affect the return “label”.
- symm_orblist of 2d array
the symmetry adapted basis
- mo2d array
the orbitals to label
- Returns:
list of symbols or integers to represent the irreps for the given orbitals
Examples:
>>> from pyscf import gto, scf, symm >>> mol = gto.M(atom='H 0 0 0; H 0 0 1', basis='ccpvdz',verbose=0, symmetry=1) >>> mf = scf.RHF(mol) >>> mf.kernel() >>> symm.label_orb_symm(mol, mol.irrep_name, mol.symm_orb, mf.mo_coeff) ['Ag', 'B1u', 'Ag', 'B1u', 'B2u', 'B3u', 'Ag', 'B2g', 'B3g', 'B1u'] >>> symm.label_orb_symm(mol, mol.irrep_id, mol.symm_orb, mf.mo_coeff) [0, 5, 0, 5, 6, 7, 0, 2, 3, 5]
- pyscf.symm.addons.route(target, nelec, orbsym)[source]#
Pick orbitals to form a determinant which has the right symmetry. If solution is not found, return []
- pyscf.symm.addons.std_symb(gpname)[source]#
std_symb(‘d2h’) returns D2h; std_symb(‘D2H’) returns D2h
- pyscf.symm.addons.symmetrize_multidim(mol, mo, s=None, check=True, tol=1e-10, keep_phase=True)[source]#
Symmetrize orbitals with respect to multidimensional irreps.
Make coefficients of partner functions of multidimensional irreps to be the same. The functions uses the convention of the libmsym interface, that introduces underscores to the labels of multidimensional irreps partners.
- Args:
- molan instance of
Mole
Symmetry-adapted basis with multidimensional irreps should be generated by libmsym
- mo2D float array
The orbital space to symmetrize
- molan instance of
- Kwargs:
- s2D float array
Overlap matrix. If not given, overlap is computed with the input mol.
- checkbool
Whether to check orthogonality of input orbitals and try to fix it
- tolfloat
Orthogonality tolerance
- keep_phasebool
Whether to keep original orbital phases, rather then make them coherent with the first partner
- Returns:
2D orbital coefficients
- pyscf.symm.addons.symmetrize_orb(mol, mo, orbsym=None, s=None, check=False)[source]#
Symmetrize the given orbitals.
This function is different to the
symmetrize_space()
: In this function, each orbital is symmetrized by removing non-symmetric components.symmetrize_space()
symmetrizes the entire space by mixing different orbitals.Note this function might return non-orthogonal orbitals. Call
symmetrize_space()
to find the symmetrized orbitals that are close to the given orbitals.- Args:
- mo2D float array
The orbital space to symmetrize
- Kwargs:
- orbsyminteger list
Irrep id for each orbital. If not given, the irreps are guessed by calling
label_orb_symm()
.- s2D float array
Overlap matrix. If given, use this overlap than the the overlap of the input mol.
- Returns:
2D orbital coefficients
Examples:
>>> from pyscf import gto, symm, scf >>> mol = gto.M(atom = 'C 0 0 0; H 1 1 1; H -1 -1 1; H 1 -1 -1; H -1 1 -1', ... basis = 'sto3g') >>> mf = scf.RHF(mol).run() >>> mol.build(0, 0, symmetry='D2') >>> mo = symm.symmetrize_orb(mol, mf.mo_coeff) >>> print(symm.label_orb_symm(mol, mol.irrep_name, mol.symm_orb, mo)) ['A', 'A', 'B1', 'B2', 'B3', 'A', 'B1', 'B2', 'B3']
- pyscf.symm.addons.symmetrize_space(mol, mo, s=None, check=True, tol=1e-09, clean=False)[source]#
Symmetrize the given orbital space.
This function is different to the
symmetrize_orb()
: In this function, the given orbitals are mixed to reveal the symmetry;symmetrize_orb()
projects out non-symmetric components for each orbital.- Args:
mol : an instance of
Mole
- mo2D float array
The orbital space to symmetrize
- Kwargs:
- s2D float array
Overlap matrix. If not given, overlap is computed with the input mol.
- checkbool
Whether to check orthogonality of input orbitals and try to fix it
- tolfloat
Orthogonality tolerance
- cleanbool
Whether to zero out symmetry forbidden orbital coefficients
- Returns:
2D orbital coefficients
Examples:
>>> from pyscf import gto, symm, scf >>> mol = gto.M(atom = 'C 0 0 0; H 1 1 1; H -1 -1 1; H 1 -1 -1; H -1 1 -1', ... basis = 'sto3g') >>> mf = scf.RHF(mol).run() >>> mol.build(0, 0, symmetry='D2') >>> mo = symm.symmetrize_space(mol, mf.mo_coeff) >>> print(symm.label_orb_symm(mol, mol.irrep_name, mol.symm_orb, mo)) ['A', 'A', 'A', 'B1', 'B1', 'B2', 'B2', 'B3', 'B3']
pyscf.symm.basis module#
Generate symmetry adapted basis
pyscf.symm.cg module#
pyscf.symm.geom module#
- exception pyscf.symm.geom.RotationAxisNotFound[source]#
Bases:
PointGroupSymmetryError
- class pyscf.symm.geom.SymmSys(atoms, basis=None)[source]#
Bases:
object
- property atom_coords#
- pyscf.symm.geom.alias_axes(axes, ref)[source]#
Rename axes, make it as close as possible to the ref axes
- pyscf.symm.geom.check_given_symm(gpname, atoms, basis=None)#
Check whether the declared symmetry (gpname) exists in the system
If basis is specified, this function checks also the basis functions have the required symmetry.
- Args:
- gpname: str
point group name
- atoms: list
[[symbol, [x, y, z]], [symbol, [x, y, z]], …]
- pyscf.symm.geom.check_symm(gpname, atoms, basis=None)[source]#
Check whether the declared symmetry (gpname) exists in the system
If basis is specified, this function checks also the basis functions have the required symmetry.
- Args:
- gpname: str
point group name
- atoms: list
[[symbol, [x, y, z]], [symbol, [x, y, z]], …]
- pyscf.symm.geom.detect_symm(atoms, basis=None, verbose=2)[source]#
Detect the point group symmetry for given molecule.
Return group name, charge center, and nex_axis (three rows for x,y,z)
- pyscf.symm.geom.rotation_mat(vec, theta)[source]#
rotate angle theta along vec new(x,y,z) = R * old(x,y,z)
- pyscf.symm.geom.subgroup(gpname, axes)#
pyscf.symm.msym module#
pyscf.symm.param module#
- pyscf.symm.param.SPHERIC_GTO_PARITY_ODD = (((0, 0, 0),), ((1, 0, 0), (0, 1, 0), (0, 0, 1)), ((1, 1, 0), (0, 1, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0)), ((0, 1, 0), (1, 1, 1), (0, 1, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0)), ((1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0)), ((0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0)), ((1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0)), ((0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0)), ((1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0)), ((0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0)), ((1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0)), ((0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0)), ((1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0)), ((0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0)), ((1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (1, 1, 0), (0, 1, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0), (1, 0, 1), (0, 0, 0)), ((0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (1, 1, 1), (0, 1, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0), (0, 0, 1), (1, 0, 0)))#
- def parity_lm(c, l):
n = 0 for lx in reversed(range(0, l+1)):
- for ly in reversed(range(0, l+1 - lx)):
lz = l - lx - ly if c[n] != 0:
return (lx%2, ly%2, lz%2)
n += 1
[[parity_lm(c, l) for c in gto.cart2sph(l).T] for l in range(16)]
pyscf.symm.sph module#
Spherical harmonics
- pyscf.symm.sph.multipoles(r, lmax, reorder_dipole=True)[source]#
Compute all multipoles upto lmax
rad = numpy.linalg.norm(r, axis=1) ylms = real_ylm(r/rad.reshape(-1,1), lmax) pol = [rad**l*y for l, y in enumerate(ylms)]
- Kwargs:
- reorder_pbool
sort dipole to the order (x,y,z)
- pyscf.symm.sph.real2spinor(l, reorder_p=True)#
- pyscf.symm.sph.real2spinor_whole(mol)#
Transformation matrix that transforms real-spherical GTOs to spinor GTOs for all basis functions
Examples:
>>> from pyscf import gto >>> from pyscf.symm import sph >>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='ccpvtz') >>> ca, cb = sph.sph2spinor_coeff(mol) >>> s0 = mol.intor('int1e_ovlp_spinor') >>> s1 = ca.conj().T.dot(mol.intor('int1e_ovlp_sph')).dot(ca) >>> s1+= cb.conj().T.dot(mol.intor('int1e_ovlp_sph')).dot(cb) >>> print(abs(s1-s0).max()) >>> 6.66133814775e-16
- pyscf.symm.sph.real_sph_vec(r, lmax, reorder_p=False)[source]#
Computes (all) real spherical harmonics up to the angular momentum lmax
- pyscf.symm.sph.sph2spinor_coeff(mol)[source]#
Transformation matrix that transforms real-spherical GTOs to spinor GTOs for all basis functions
Examples:
>>> from pyscf import gto >>> from pyscf.symm import sph >>> mol = gto.M(atom='H 0 0 0; F 0 0 1', basis='ccpvtz') >>> ca, cb = sph.sph2spinor_coeff(mol) >>> s0 = mol.intor('int1e_ovlp_spinor') >>> s1 = ca.conj().T.dot(mol.intor('int1e_ovlp_sph')).dot(ca) >>> s1+= cb.conj().T.dot(mol.intor('int1e_ovlp_sph')).dot(cb) >>> print(abs(s1-s0).max()) >>> 6.66133814775e-16
- pyscf.symm.sph.sph_pure2real(l, reorder_p=True)[source]#
Transformation matrix: from the pure spherical harmonic functions Y_m to the real spherical harmonic functions O_m:
O_m = \sum Y_m' * U(m',m)
Y(-1) = 1/sqrt(2){-iO(-1) + O(1)}; Y(1) = 1/sqrt(2){-iO(-1) - O(1)} Y(-2) = 1/sqrt(2){-iO(-2) + O(2)}; Y(2) = 1/sqrt(2){iO(-2) + O(2)} O(-1) = i/sqrt(2){Y(-1) + Y(1)}; O(1) = 1/sqrt(2){Y(-1) - Y(1)} O(-2) = i/sqrt(2){Y(-2) - Y(2)}; O(2) = 1/sqrt(2){Y(-2) + Y(2)}
- Kwargs:
reorder_p (bool): Whether the p functions are in the (x,y,z) order.
- Returns:
2D array U_{complex,real}
Module contents#
This module offers the functions to detect point group symmetry, basis symmetrization, Clebsch-Gordon coefficients. This module works as a plugin of PySCF package. Symmetry is not hard coded in each method.