qsymm.groups module

class qsymm.groups.ContinuousGroupGenerator(R=None, U=None)[source]

Bases: object

A generator of a continuous group.

Generates a family of symmetry operators that act on the Hamiltonian as:

\[H(k) → \exp{-iλU} H(\exp{iλR} k) \exp{iλU}\]

with λ a real parameter.

Parameters:
  • R (ndarray (optional)) – Real space rotation generator: Hermitian antisymmetric. If not provided, the zero matrix is used.

  • U (ndarray (optional)) – Hilbert space unitary rotation generator: Hermitian. If not provided, the zero matrix is used.

R
U
apply(model)[source]

Return copy of model H(k) with applied infinitesimal generator. 1j * (H(k) U - U H(k)) + 1j * dH(k)/dk_i R_{ij} k_j

qsymm.groups.L_matrices(d=3, l=1)[source]

Construct real space rotation generator matrices in d=2 or 3 dimensions. Can also be used to get angular momentum operators for real atomic orbitals in 3 dimensions, for p-orbitals use l=1, for d-orbitals l=2. The basis of p-orbitals is p_x, p_y, p_z, for d-orbitals d_{x^2 - y^2}, d_{3 z^2 - r^2}, d_{xy}, d_{yz}, d_{zx}. The matrices are all purely imaginary and antisymmetric. To generate finite rotations, use ‘spin_rotation(n, L)’.

class qsymm.groups.PointGroupElement(R, conjugate=False, antisymmetry=False, U=None, _strict_eq=False, *, locals=None)[source]

Bases: object

An element of a point group.

Parameters:
  • R (sympy.ImmutableMatrix or array) – Real space rotation action of the operator. Square matrix with size of the number of spatial dimensions.

  • conjugate (boolean (default False)) – Whether the operation includes conplex conjugation (antiunitary operator)

  • antisymmetry (boolean (default False)) – Whether the operator flips the sign of the Hamiltonian (antisymmetry)

  • U (array, str, SymPy expression, or None (default)) – The unitary action on the Hilbert space. May be None, to be able to treat symmetry candidates

  • _strict_eq (boolean (default False)) – Whether to test the equality of the unitary parts when comparing with other PointGroupElements. By default the unitary parts are ignored. If True, PointGroupElements are considered equal, if the unitary parts are proportional, an overall phase difference is still allowed.

  • locals (dict or None (default)) – Additional namespace entries for sympify. May be used to simplify input of matrices or modify input before proceeding further. For example: locals={'sigma_plus': [[0, 2], [0, 0]]}.

Notes

As U is floating point and has a phase ambiguity at least, hence it is ignored when comparing objects by default.

R is the real space rotation acion. Do not include minus sign for the k-space action of antiunitary operators, such as time reversal. This minus sign will be included automatically if ‘conjugate=True’.

For most uses R can be provided as a floating point array. It is necessary to use exact sympy matrix representation if the PGE has to act on Models with complicated momentum dependence (not polynomial), as the function parts of models are compared exactly. If the momentum dependence is periodic (sine, cosine and exponential), use BlochModel, this works with floating point rotations.

R
U
antisymmetry
apply(model)[source]

Return copy of model with applied symmetry operation.

if unitary: (+/-) U H(inv(R) k) U^dagger if antiunitary: (+/-) U H(- inv(R) k).conj() U^dagger

(+/-) stands for (symmetry / antisymmetry)

If self.U is None, U is taken as the identity.

conjugate
identity()[source]

Return identity element with the same structure as self.

inv()[source]

Invert PointGroupElement

class qsymm.groups.PrettyList(iterable=(), /)[source]

Bases: list

Subclass of list that displays its elements with latex.

qsymm.groups.chiral(realspace_dim, U=None)[source]

Return a chiral symmetry operator

Parameters:
  • realspace_dim (int) – Realspace dimension

  • U (ndarray (optional)) – The unitary action on the Hilbert space. May be None, to be able to treat symmetry candidates

Returns:

P

Return type:

PointGroupElement

qsymm.groups.cubic(tr=True, ph=True, generators=False, spin=None)[source]

Generate cubic point group in standard basis.

Parameters:
  • tr (bool (default True)) – Whether to include time-reversal and particle-hole symmetry.

  • ph (bool (default True)) – Whether to include time-reversal and particle-hole symmetry.

  • generators (bool (default false)) – Only return the group generators if True.

  • spin (float or sequence of arrays (optional)) – Spin representation to use for the unitary action of the operator. If not provided, the PointGroupElements have the unitary action set to None. If float is provided, it should be integer or half-integer specifying the spin representation in the standard basis, see spin_matrices. Otherwise a sequence of 3 arrays of identical square size must be provided representing 3 components of the angular momentum operator. The unitary action of rotation operator is U = exp(-i n⋅s). If ph is True, spin may not be provided, as it is not possible to deduce the unitary representation of particle-hole symmetry from spin alone. In this case construct the particle-hole operator manually.

Returns:

  • set of PointGroupElement objects with integer rotations

  • Notes

  • —— – We assume inversion acts trivially in spin space.

qsymm.groups.generate_group(gens)[source]

Generate group from gens

Parameters:

gens (iterable of PointGroupElement) – generator set of the group

Returns:

group – group generated by gens, closed under multiplication

Return type:

set of PointGroupElement

qsymm.groups.generate_subgroups(group)[source]

Generate all subgroups of group, including the trivial group and itself.

Parameters:

group (set of PointGroupElement) – A closed group as set of its elements.

Returns:

subgroups – frozesets are the subgroups, sets are a generator set of the subgroup.

Return type:

dict of forzenset: set

qsymm.groups.hexagonal(dim=2, tr=True, ph=True, generators=False, sympy_R=True, spin=None)[source]

Generate hexagonal point group in standard basis in 2 or 3 dimensions. Mirror symmetries with the main coordinate axes as normals are included. In 3D the hexagonal axis is the z axis.

Parameters:
  • dim (int (default 2)) – Real sapce dimensionality, 2 or 3.

  • tr (bool (default True)) – Whether to include time-reversal and particle-hole symmetry.

  • ph (bool (default True)) – Whether to include time-reversal and particle-hole symmetry.

  • generators (bool (default True)) – Only return the group generators if True.

  • sympy_R (bool (default True)) – Whether the rotation matrices should be exact sympy representations.

  • spin (float or sequence of arrays (optional)) – Spin representation to use for the unitary action of the operator. If not provided, the PointGroupElements have the unitary action set to None. If float is provided, it should be integer or half-integer specifying the spin representation in the standard basis, see spin_matrices. Otherwise a sequence of 3 arrays of identical square size must be provided representing 3 components of the angular momentum operator. The unitary action of rotation operator is U = exp(-i n⋅s). In 2D the z axis is assumed to be the rotation axis. If ph is True, spin may not be provided, as it is not possible to deduce the unitary representation of particle-hole symmetry from spin alone. In this case construct the particle-hole operator manually.

Returns:

  • set of PointGroupElements

  • Notes

  • —— – Warning: in 2D the real space action of a mirror and and a 2-fold rotation around an axis in the plane is identical, however the action on angular momentum is different. Here we consider the action of the mirror, which is the same as the action of a 2-fold rotation around the mirror axis, assuming inversion acts trivially.

qsymm.groups.identity(dim, shape=None)[source]

Return identity operator with appropriate shape.

Parameters:
  • dim (int) – Dimension of real space.

  • shape (int (optional)) – Size of the unitary part of the operator. If not provided, U is set to None.

Returns:

id

Return type:

PointGroupElement

qsymm.groups.inversion(realspace_dim, U=None)[source]

Return an inversion operator

Parameters:
  • realspace_dim (int) – Realspace dimension

  • U (ndarray (optional)) – The unitary action on the Hilbert space. May be None, to be able to treat symmetry candidates

Returns:

P

Return type:

PointGroupElement

qsymm.groups.is_sympy_matrix(R)[source]
qsymm.groups.mirror(axis, U=None, spin=None)[source]

Return a mirror operator

Parameters:
  • axis (ndarray) – Normal of the mirror. The dimensionality of the operator is the same as the length of axis.

  • U (ndarray (optional)) – The unitary action on the Hilbert space. May be None, to be able to treat symmetry candidates

  • spin (float or sequence of arrays (optional)) – Spin representation to use for the unitary action of the operator. If float is provided, it should be integer or half-integer specifying the spin representation in the standard basis, see spin_matrices. Otherwise a sequence of 3 arrays of identical square size must be provided representing 3 components of the angular momentum operator. The unitary action of mirror operator is U = exp(-i π n⋅s) where n is normalized to 1. In 2D the axis is treated as x and y coordinates. Only one of U and spin may be provided.

Returns:

  • P (PointGroupElement)

  • Notes

  • —— – Warning: in 2D the real space action of a mirror and and a 2-fold rotation around an axis in the plane is identical, however the action on angular momentum is different. Here we consider the action of the mirror, which is the same as the action of a 2-fold rotation around the mirror axis.

qsymm.groups.particle_hole(realspace_dim, U=None)[source]

Return a particle-hole symmetry operator

Parameters:
  • realspace_dim (int) – Realspace dimension

  • U (ndarray (optional)) – The unitary action on the Hilbert space. May be None, to be able to treat symmetry candidates

Returns:

P

Return type:

PointGroupElement

qsymm.groups.pretty_print_cgg(g, latex=False)[source]

Return a human readable string representation of ContinuousGroupGenerator

Parameters:
  • g (ContinuousGroupGenerator) – Point group element to be represented.

  • latex (bool (default False)) – Whether to output LateX formatted string.

Returns:

name – The first line is the action on the Hamiltonian, the following lines display the real space rotation as R(ϕ) for 2D rotation and R(ϕ, axis) for 3D rotation (axis is not normalized), and the conserved quantity as a matrix. If either of these is trivial, it is omitted. Note that a ContinuousGroupGenerator defines a continuous group of symmetries, so there is always a free parameter ϕ.

Return type:

string

qsymm.groups.pretty_print_pge(g, full=False, latex=False)[source]

Return a human readable string representation of PointGroupElement

Parameters:
  • g (PointGroupElement) – Point group element to be represented.

  • full (bool (default False)) – Whether to return a full representation. The default short representation only contains the real space action and the symbol of the Altland-Zirnbauer part of the symmetry (see below). The full representation presents the symmetry action on the Hamiltonian and the unitary Hilbert-space action if set.

  • latex (bool (default False)) – Whether to output LateX formatted string.

Returns:

name – In the short representation it is a sting rot_name + az_name. In the long representation the first line is the action on the Hamiltonian, the second line is rot_name and the third line is the unitary action as a matrix, if set.

rot_name can be: - 1 for identity - I for inversion (in 1D mirror is the same as inversion) - R(angle) for 2D rotation - R(angle, axis) for 3D rotation (axis is not normalized) - M(normal) for mirror - S(angle, axis) for 3D rotoinversion (axis is not normalized)

az_name can be: - T for time-reversal (antiunitary symmetry) - P for particle-hole (antiunitary antisymmetry) - C for chiral (unitary antisymmetry) - missing if the symmetry is unitary

Return type:

string

qsymm.groups.rotation(angle, axis=None, inversion=False, U=None, spin=None)[source]

Return a rotation operator

Parameters:
  • angle (float) – Rotation angle in units of 2 pi.

  • axis (ndarray or None (default)) – Rotation axis, optional. If not provided, a 2D rotation is generated around the axis normal to the plane. If a 3D vector is provided, a 3D rotation is generated around this axis. Does not need to be normalized to 1.

  • inversion (bool (default False)) – Whether to generate a rotoinversion. By default a proper rotation is returned. Only valid in 3D.

  • U (ndarray (optional)) – The unitary action on the Hilbert space. May be None, to be able to treat symmetry candidates

  • spin (float or sequence of arrays (optional)) – Spin representation to use for the unitary action of the operator. If float is provided, it should be integer or half-integer specifying the spin representation in the standard basis, see spin_matrices. Otherwise a sequence of 3 arrays of identical square size must be provided representing 3 components of the angular momentum operator. The unitary action of rotation operator is U = exp(-i n⋅s). In 2D the z axis is assumed to be the rotation axis. Only one of U and spin may be provided.

Returns:

P

Return type:

PointGroupElement

qsymm.groups.rotation_to_angle(R)[source]
qsymm.groups.set_multiply(G, H)[source]
qsymm.groups.spin_matrices(s, include_0=False)[source]

Construct spin-s matrices for any half-integer spin.

Parameters:
  • s (float or int) – Spin representation to use, must be integer or half-integer.

  • include_0 (bool (default False)) – If include_0 is True, S[0] is the identity, indices 1, 2, 3 correspond to x, y, z. Otherwise indices 0, 1, 2 are x, y, z.

Returns:

Sequence of spin-s operators in the standard spin-z basis. Array of shape (3, 2*s + 1, 2*s + 1), or if include_0 is True (4, 2*s + 1, 2*s + 1).

Return type:

ndarray

qsymm.groups.spin_rotation(n, s, roundint=False)[source]

Construct the unitary spin rotation matrix for rotation specified by the vector n (in radian units) with angular momentum s, given by U = exp(-i n⋅s).

Parameters:
  • n (iterable) – Rotation vector. Its norm is the rotation angle in radians.

  • s (float or sequence of arrays) – Spin representation to use for the unitary action of the operator. If float is provided, it should be integer or half-integer specifying the spin representation in the standard basis, see spin_matrices. Otherwise a sequence of 3 arrays of identical square size must be provided representing 3 components of the angular momentum operator.

  • roundint (bool (default False)) – If roundint is True, result is converted to integer tinyarray if possible.

Returns:

U – Unitary spin rotation matrix of the same shape as the spin matrices or (2*s + 1, 2*s + 1).

Return type:

ndarray

qsymm.groups.square(tr=True, ph=True, generators=False, spin=None)[source]

Generate square point group in standard basis.

Parameters:
  • tr (bool (default True)) – Whether to include time-reversal and particle-hole symmetry.

  • ph (bool (default True)) – Whether to include time-reversal and particle-hole symmetry.

  • generators (bool (default false)) – Only return the group generators if True.

  • spin (float or sequence of arrays (optional)) – Spin representation to use for the unitary action of the operator. If not provided, the PointGroupElements have the unitary action set to None. If float is provided, it should be integer or half-integer specifying the spin representation in the standard basis, see spin_matrices. Otherwise a sequence of 3 arrays of identical square size must be provided representing 3 components of the angular momentum operator. The unitary action of rotation operator is U = exp(-i n⋅s). In 2D the z axis is assumed to be the rotation axis. If ph is True, spin may not be provided, as it is not possible to deduce the unitary representation of particle-hole symmetry from spin alone. In this case construct the particle-hole operator manually.

Returns:

  • set of PointGroupElement objects with integer rotations

  • Notes

  • —— – Warning: in 2D the real space action of a mirror and and a 2-fold rotation around an axis in the plane is identical, however the action on angular momentum is different. Here we consider the action of the mirror, which is the same as the action of a 2-fold rotation around the mirror axis, assuming inversion acts trivially.

qsymm.groups.symmetry_from_permutation(R, perm, norbs, onsite_action=None, antiunitary=False, antisymmetry=False)[source]

Construct symmetry operator for lattice systems with multiple sites.

Parameters:
  • R (real space rotation) –

  • perm (dict : {site: image_site}) – permutation of the sites under the symmetry action

  • norbs (OrderedDict : {site : norbs_site} or tuple of tuples ((site, norbs_site), )) – sites are ordered in the order specified, with blocks of size norbs_site corresponding to each site.

  • onsite_action (dict : {site: ndarray} or ndarray or None) – onsite symmetry action, such as spin rotation for each site. If only one array is specified, it is used for every site. If None (default), identity is used on every site. Size of the arrays must be consistent with norbs.

  • antiunitary (bool) –

  • antisymmetry (bool) –

Returns:

  • g (PointGroupElement) – PointGroupElement corresponding to the operation.

  • Notes

  • ——

  • Sites can be indexed by any hashable identifiers, such as integers or stings.

qsymm.groups.time_reversal(realspace_dim, U=None, spin=None)[source]

Return a time-reversal symmetry operator

Parameters:
  • realspace_dim (int) – Realspace dimension

  • U (ndarray (optional)) – The unitary action on the Hilbert space. May be None, to be able to treat symmetry candidates.

  • spin (float or sequence of arrays (optional)) – Spin representation to use for the unitary action of the time reversal operator. If float is provided, it should be integer or half-integer specifying the spin representation in the standard basis, see spin_matrices. Otherwise a sequence of 3 arrays of identical square size must be provided representing 3 components of the angular momentum operator. The unitary action of time-reversal operator is U = exp(-i π s_y). Only one of U and spin may be provided.

Returns:

T

Return type:

PointGroupElement