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)

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)

Wigner small-d matrix (in z-y-z convention)

pyscf.symm.Dmatrix.get_euler_angles(c1, c2)

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.direct_prod(orbsym1, orbsym2, groupname='D2h')
pyscf.symm.addons.eigh(h, orbsym)

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.irrep_id2name(gpname, irrep_id)

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 sybmol, str

pyscf.symm.addons.irrep_name(pgname, irrep_id)
pyscf.symm.addons.irrep_name2id(gpname, symb)

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)

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)

Pick orbitals to form a determinant which has the right symmetry. If solution is not found, return []

pyscf.symm.addons.std_symb(gpname)

std_symb(‘d2h’) returns D2h; std_symb(‘D2H’) returns D2h

pyscf.symm.addons.symmetrize_orb(mol, mo, orbsym=None, s=None, check=False)

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-orthorgonal 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-07)

Symmetrize the given orbital space.

This function is different to the symmetrize_orb(): In this function, the given orbitals are mixed to reveal the symmtery; symmetrize_orb() projects out non-symmetric components for each orbital.

Args:
mo2D float array

The orbital space to symmetrize

Kwargs:
s2D float array

Overlap matrix. If not given, overlap is computed with 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_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.basis.dump_symm_adapted_basis(mol, so)
pyscf.symm.basis.linearmole_irrep_id2symb(gpname, irrep_id)
pyscf.symm.basis.linearmole_irrep_symb2id(gpname, symb)
pyscf.symm.basis.linearmole_symm_adapted_basis(mol, gpname, orig=0, coordinates=None)
pyscf.symm.basis.linearmole_symm_descent(gpname, irrepid)

Map irreps to D2h or C2v

pyscf.symm.basis.symm_adapted_basis(mol, gpname, orig=0, coordinates=None)
pyscf.symm.basis.symmetrize_matrix(mat, so)
pyscf.symm.basis.tot_parity_odd(op, l, m)

pyscf.symm.cg module

pyscf.symm.cg.cg_spin(l, jdouble, mjdouble, spin)

Clebsch Gordon coefficient of <l,m,1/2,spin|j,mj>

pyscf.symm.geom module

exception pyscf.symm.geom.RotationAxisNotFound

Bases: RuntimeError

class pyscf.symm.geom.SymmSys(atoms, basis=None)

Bases: object

cartesian_tensor(n)
has_icenter()
has_improper_rotation(axis, n)
has_mirror(perp_vec)
has_rotation(axis, n)
search_c2x(zaxis, n)

C2 axis which is perpendicular to z-axis

search_c_highest(zaxis=None)
search_mirrorx(zaxis, n)

mirror which is parallel to z-axis

search_possible_rotations(zaxis=None)

If zaxis is given, the rotation axis is parallel to zaxis

symmetric_for(op)
pyscf.symm.geom.alias_axes(axes, ref)

Rename axes, make it as close as possible to the ref axes

pyscf.symm.geom.argsort_coords(coords, decimals=None)
pyscf.symm.geom.as_subgroup(topgroup, axes, subgroup=None)
pyscf.symm.geom.check_given_symm(gpname, atoms, basis=None)
pyscf.symm.geom.closest_axes(axes, ref)
pyscf.symm.geom.detect_symm(atoms, basis=None, verbose=2)

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.get_subgroup(gpname, axes)
pyscf.symm.geom.householder(vec)
pyscf.symm.geom.parallel_vectors(v1, v2, tol=1e-05)
pyscf.symm.geom.rotation_mat(vec, theta)

rotate angle theta along vec new(x,y,z) = R * old(x,y,z)

pyscf.symm.geom.shift_atom(atoms, orig, axis)
pyscf.symm.geom.sort_coords(coords, decimals=None)
pyscf.symm.geom.subgroup(gpname, axes)
pyscf.symm.geom.symm_identical_atoms(gpname, atoms)

Requires

pyscf.symm.geom.symm_ops(gpname, axes=None)

pyscf.symm.param module

pyscf.symm.sph module

Spherical harmonics

pyscf.symm.sph.cart2spinor(l)

Cartesian to spinor for angular moment l

pyscf.symm.sph.multipoles(r, lmax, reorder_dipole=True)

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)

Computes (all) real spherical harmonics up to the angular momentum lmax

pyscf.symm.sph.sph2spinor(l, reorder_p=True)
pyscf.symm.sph.sph2spinor_coeff(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.sph_pure2real(l, reorder_p=True)

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}

pyscf.symm.sph.sph_real2pure(l, reorder_p=True)

Transformation matrix: from real spherical harmonic functions to the pure spherical harmonic functions.

Kwargs:

reorder_p (bool): Whether the real p functions are in the (x,y,z) order.

Module contents

This module offers the functions to detect point group symmetry, basis symmetriziation, Clebsch-Gordon coefficients. This module works as a plugin of PySCF package. Symmetry is not hard coded in each method.