Package documentation#

Main module#

MHI – “Multi-Hadron Interpolators”

Module for constructing block-diagonalization / change-of-basis matrices to map products of N local plane-wave operators into irreps of the cubic group. Includes appropriate generalizations for spin and internal symmetries.

mhi.mhi.mhi(momenta, spin_irreps=None, internal_symmetry=None, verbose=False)[source]#

General-purpose driver function for construction of change-of-basis / block-diagonalization matrices which project linear combinations of plane- wave states onto irreps of the cubic group.

Parameters:
momenta(nparticles, 3) or (3, ) array_like

The ordered momenta.

particle_namesarray_like

The particle names, e.g., [‘n’, ‘p’].

verbosebool

Whether or not to print extra diagnostic information.

return_Dmmbool

Whether or not to return the momentum-(spin) representation matrices.

internal_symmetrylist of WeightedPermutation

The exchange group projector defined in the group algebra.

Returns:
resultIrrepDecomposition

object containing the results of the irrep decomposition.

Notes

The algorithm is as follows:
  • Compute the little group of the total momentum

  • Compute irrep matrices of the little group, \(D_{\mu\nu}(R)\)

  • Compute the momentum(-spin) representation matrices, \(D_{mm'}(R)\)

  • Apply exchange-group projection, giving \(\hat{D}(R) = P D(R) P\)

  • Apply Schur’s algorithm and transition operators to construct the block-diagonalization matrices

class mhi.mhi.WeightedPermutation(weight, perm)[source]#

A complex scalar weight times a permutation group element.

Attributes:
weightcomplex or float

Scalar weight multiplying the permutation.

perm(n,) ndarray

Permutation expressed as an array.

class mhi.mhi.Isomorphism(g, perm)[source]#

The group element and permutation specifying a subgroup isomorphism.

Attributes:
g(n, n) ndarray

Group element used to conjugate the subgroup as \(g \cdot H \cdot g^{-1}\)

perm(|H|,) ndarray

Permutation p mapping from conjugated elements to target subgroup as \(H' = g \cdot H \cdot g^{-1}[p]\)

class mhi.mhi.SpinShellTuple(momenta, spins)[source]#

Pairing of momenta and spins making up an orbit with non-zero spins.

Attributes:
momenta(norbit, nmomenta, 3) ndarray

List of momentum lists in the orbit.

spins(norbit, nspin) ndarray

List of spin configurations in the orbit.

mhi.mhi.levi_civita(dim)[source]#

The totally antisymmetric Levi-Civita tensor in arbitrary dimensions.

Parameters:
dimint

The number of dimensions.

Returns:
arrndarray

The Levi-Civita tensor.

Notes

Implementation from StackOverflow user JGibs: https://stackoverflow.com/questions/59740966/levi-civita-tensor-in-numpy

mhi.mhi.unique_permutations(seq)[source]#

Yield only unique permutations of seq in an efficient way.

Parameters:
seq: array_like

The elements to be permuted

Returns:
seq: array_like

The permuated sequences

Notes

A python implementation of Knuth’s “Algorithm L”, also known from the std::next_permutation function of C++, and as the permutation algorithm of Narayana Pandita.

Code taken from a post by StackOverflow user Lauritz V. Thaulow: https://stackoverflow.com/questions/12836385/how-can-i-interleave-or-create-unique-permutations-of-two-strings-without-recur/12837695#12837695

Examples

>>> for perm in unique_permutations([1, 1, 2]):
>>>     print(perm)
[1, 1, 2]
[1, 2, 1]
[2, 1, 1]
>>> for perm in unique_permutations(['a', 'a', 'a', 'b', 'b']):
>>>     print(perm)
['a', 'a', 'a', 'b', 'b']
['a', 'a', 'b', 'a', 'b']
['a', 'a', 'b', 'b', 'a']
['a', 'b', 'a', 'a', 'b']
['a', 'b', 'a', 'b', 'a']
['a', 'b', 'b', 'a', 'a']
['b', 'a', 'a', 'a', 'b']
['b', 'a', 'a', 'b', 'a']
['b', 'a', 'b', 'a', 'a']
['b', 'b', 'a', 'a', 'a']
mhi.mhi.symmetrize(arr)[source]#

Computes the completely symmetrized version of the input array.

Parameters:
arrarray_like

A tensor of arbitrary rank.

Returns:
arr_symndarray

The completely symmetrized version of the input tensor.

Notes

In index notation symmetrization is the map which sends \(T_{ij...n} \to T_{(ij...n)}\), where parantheses denote symmetrization. For instance, the symmetrization of a three-index tensor \(T_{ijk}\) is \(T_{(ijk)} = \tfrac{1}{6}(T_{ijk} + T_{ikj} + T_{jik} + T_{jki} + T_{kij} + T_{kji})\)

mhi.mhi.tensor_product(a, b)[source]#

Computes the tensor product between tensors a and b.

Parameters:
a, barray_like
Returns:
tensorarray_like

The tensor product of arrays a and b.

Notes

In index notation, this fuction computes the output tensor T defined by \(T_{ij \dots k rs \dots t} = a_{ij \dots k} b_{rs \dots t}\).

mhi.mhi.tensor_nfold(*tensors)[source]#

Computes the n-fold tensor product between all the tensors in a list.

Parameters:
tensors(variable number of) ndarray

The input tensors for the product.

Returns:
tensorndarray

The product of all the tensors.

Notes

Suppose the input tensors are \(\{a_i, b_{jk}, c_{lmn}\}\). This function computes the output tensor T defined by \(T_{ijklmn} = a_i b_{jk} c_{lmn}\).

mhi.mhi.decompose(monomial)[source]#

Decomposes a sympy monomial of the form c * x**nx * y**ny * z**nz into a coefficient c and a triplet of exponents (nx, ny, nz).

Parameters:
monomialsympy.core.mul.Mul

A monomial in the variables {‘x’, ‘y’, ‘z’}

Returns:
(coeff, exponents)(complex, tuple)

The coefficient and exponents (nx, ny, nz) of the monomial

mhi.mhi.compute_polarization(monomial)[source]#

Computes the polarization tensor of a given monomial.

Parameters:
monomialsympy.core.mul.Mul

A monomial in the variables {‘x’, ‘y’, ‘z’}

Returns:
polarization_tensorndarray

The polarization, a totally symmetric tensor

mhi.mhi.contract_across(arr, vec)[source]#

Conctract a vector full across all the free indices, yielding a scalar.

Parameters:
arrarray_like

Tensor of generic rank

vecarray_like

The vector to contract against all the free indices

Returns:
scalarfloat or complex

The product c = vec[i]*vec[j]*…*vec[n]*arr[i,j,…,n]

mhi.mhi.compute_restitution(polarization)[source]#

Computes the monomial restitution of a polarization tensor.

Parameters:
polarizationarray_like

The (completely symmetric) polarization tensor

Returns:
restitutionmonomial (sympy.core.mul.Mul)

The monomial resulting from contracting the vector [x,y,z] across all the free indices of the polarization tensor.

mhi.mhi.polarize(fcn)[source]#

Computes the polarization tensor associated with a basis function.

Parameters:
fcnsympy basis function

The basis function, a sum of monomials in {‘x’, ‘y’, ‘z’}.

Returns:
arrndarray

The polarization tensor of rank d, where d is the degree of the monomials within the basis function.

mhi.mhi.transform(arr, group_element, verbose=False)[source]#

Computes the “rotation” transformation of a tensor of arbitrary rank, A[a,b,…,c] –> R[a,x] R[b,y] … R[c,z] A[x,y,…,z].

Parameters:
arr(M, M, ..., M) array_like

The tensor to transform.

group_element(M, M) array_like

The group element applying the transformation.

verbosebool

Whether to print the indices used in Einstein summation notation.

Returns:
arr_transformedndarray

The transformed tensor.

mhi.mhi.tensor_inner(a, b, verbose=False)[source]#

Computes the inner product between two tensors a[i,j,…,k]*b[i,j,…,k].

Parameters:
a(M, M, ... , M) array_like
b(M, M, ... , M) array_like
Returns:
cfloat or complex

The result of the inner product across all free indices

mhi.mhi.make_tensor_product_space(dims)[source]#

Makes a tensor product space with the specified dimensions.

Parameters:
dimslist

Integers specifying the dimensions

Returns:
tensor_product_spacendarray

Notes

These spaces arise naturally in the construction of combined “momentum- spin” orbits below. For the portion of the orbit associated with spin, it’s assumed that the spin has already been decomposed into irreps of known dimensions. By definition of the irrep space, group transformation simply permute the basis elements of the irrep space.

Examples

>>> arr = make_tensor_product_space([1,2,3])
>>> print(arr.shape)
(6, 3)
>>> print(arr)
array([
    [0, 0, 0],
    [0, 0, 1],
    [0, 0, 2],
    [0, 1, 0],
    [0, 1, 1],
    [0, 1, 2]])
class mhi.mhi.DiracPauli(verbose=False)[source]#

Wrapper for Dirac matrices in the Dirac-Pauli basis.

rotation_vec(omega)[source]#

Computes the spinor rotation matrix associated with the 3-vector omega.

Parameters:
omega: ``(3,)`` ndarray or list

The vector specifying the rotation

Returns:
arr(4, 4) ndarray

The complex roation matrix

rotation(theta, direction)[source]#

Computes the spinor rotation matrix associated with an rotation of angle “theta” around the kth axis.

Parameters:
theta: float

The rotation angle

direction: int

The axis number, with {1,2,3} <–> {x,y,z}.

Returns:
arr(4, 4) ndarray

The complex rotation matrix

mhi.mhi.get_nprod(j)[source]#

Gets the number of spin-1/2 copies present for a given half-integer j.

For example, 3/2 = 1/2 otimes 1/2 otimes 1/2 –> 3 copies

Parameters:
jint or float

The total spin.

Returns:
nint

The number of spin-1/2 copies present.

mhi.mhi.make_spinor_array(spinor)[source]#

Builds a normalized “state” with total (j, jz, parity) inside a suitable tensor product space of spin-1/2 states.

Parameters:
spinorSpinorTuple

The namedtuple specifying (j, jz, parity)

Returns:
vecndarray

The normalized vector of length 4**(2*j)

mhi.mhi.make_oh()[source]#

Constructs a presentation of the cubic group Oh with a standardized ordering of group elements.

Parameters:
None
Returns:
group(48, 3, 3) ndarray
mhi.mhi.rotation(theta, direction)[source]#

Computes the rotation matrix by the angle theta around a particular axis.

Parameters:
thetafloat

The angle

direction: int

The axis number, with {1,2,3} <–> {x,y,z}.

Returns:
arr(3, 3) ndarray

The rotation matrix

mhi.mhi.rotation_vec(omega)[source]#

Computes the rotation matrix associated with the three-vector omega.

Parameters:
omega(3, ) ndarray or list

The vector defining the rotation

Returns:
arr(3, 3) ndarray

The rotation matrix

mhi.mhi.make_ohd()[source]#

Constructs a presentation of the “spinorial” double cover OhD of the cubic group ordering of group elements.

Parameters:
None
Returns:
group(96, 4, 4) ndarray
mhi.mhi.make_spinorial_little_group(little)[source]#

Makes the “spinorial” little group associated with the double cover OhD, given a little group G in Oh.

Parameters:
little(|G|, 3, 3) array_like

The little group G.

Returns:
group(2*|G|, 4, 4) ndarray

The double-cover little group.

mhi.mhi.make_stabilizer(momenta, group)[source]#

Constructs the stabilizer subgroup of an ordered set of momenta. The stabilizer group is the subgroup that leaves the ordered set invariant.

Parameters:
momenta(nmomenta, 3) or (3, ) array_like

The ordered set of momenta that must be left invariant.

group(|G|, 3, 3) array_like

The total group.

Returns:
stabilizer(|H|, 3, 3) ndarray

The stabilizer group H.

mhi.mhi.make_canonical_stabilizer(name, group)[source]#

Computes the stabilizer subgroup with a canonical orientation inside the larger group Oh.

Parameters:
namestr

The name of the desired stabilizer subgroup

group(|G|, 3, 3) array_like

The group

Returns:
stabilizer(|H|, 3, 3) array_like

The stabilizer subgroup

mhi.mhi.identify_stabilizer(stabilizer)[source]#

Identifies the name of the stabilizer group “H” by checking its order.

Parameters:
stabilizer(|H|, 3, 3) array_like

The stabilizer group.

Returns:
namestr

The name of the stabilizer group.

mhi.mhi.make_little_and_stabilizer(momenta, group)[source]#

Computes the little group and stabilizer group of the ordered set of momenta.

Parameters:
momenta(nmomenta, 3) array_like

The momenta

group(|G|, 3, 3) array_like

The group

Returns:
groupstuple = (ndarray, ndarray)

A pair of groups arranged as (“little group”, “stabilizer group”)

Notes

A word on physics naming conventions. Let G be a group, let momenta be a set of ordered momenta, and let ktot be the total momentum (i.e., the sum of the momenta). The “little group” is the subgroup of G that leaves ktot invariant. The “stabilizer group” is the subgroup of G that leaves momenta invariant.

mhi.mhi.conjugate(g, h)[source]#

Computes the conjugate of group element h by group element g: \(g.h.g^{-1}\). Assumes that g is an orthogonal matrix so that \(g^{-1} = g^T\).

Parameters:
g(n, n) ndarray
h(n, n) ndarray
Returns:
h_conj(n, n) ndarray

The conjugated element, \(g.h.g^{-1}\)

mhi.mhi.conjugate_group(g, group)[source]#

Computes the conjugate of the group G by the group element g: \(g.G.g^{-1}\). Assumes that g is an orthogonal matrix so that \(g^{-1} = g^T\).

Parameters:
group(|G|, n, n) ndarray

The group G.

h(n, n) ndarray

The conjugating element.

Returns:
group_conj(|G|, n, n) ndarray

The conjugated group \(g.G.g^{-1}\)

mhi.mhi.find_subgroup_isomorphism(group, subgroup_h1, subgroup_h2)[source]#

Finds the isomorphism between conjugate subgroups H1 and H2.

Parameters:
group(|G|, n, n) ndarray

The group G.

subgroup_h1(|H|, n, n) ndarray

The subgroup H1.

subgroup_h2(|H|, n, n) ndarray

The subgroup H2.

Returns:
(g, perm)Isomorphism

The group element and permutation specifying the isomorphism.

mhi.mhi.apply_isomorphism(group, isomorphism)[source]#

Applies the isomorphism (g, perm) to the group G via \((g.G.g^{-1})[perm]\).

Parameters:
group(|G|, n, n) ndarray
isomorphismIsomorphism
Returns:
group_iso(|G|, n, n) ndarray

The group after applying the isomorphism, \((g.G.g^{-1})[perm]\)

mhi.mhi.force_hash(arr)[source]#

Computes a hash for an array.

class mhi.mhi.HashableArray(dset)[source]#

Minimal wrapper to make arrays hashable based on their value at initialization. Basically copied from the docs: https://docs.scipy.org/doc/numpy/user/basics.subclassing.html

class mhi.mhi.MomentumSpinOrbitElement(momenta, spins)[source]#

Wrapper for easy manipulations involving elements of momentum-spin orbits.

mhi.mhi.make_momentum_orbit(momenta, group, exchange_group=None)[source]#

Computes the orbit of an ordered set of vectors under a group action.

Parameters:
momenta(nmomenta, 3) array_like

The ordered momenta.

group(|G|, 3, 3) array_like

The group matrices.

exchange_grouplist of WeightedPermutation

The exchange group projector.

Returns:
orbitlist

The matrices corresponding to ordered momenta in the orbit.

Notes

In abstract algebra, an orbit is often considered as an unordered set. For numerical work, it is useful to work with orbits ordered in a standard way. The convention in the present work is that the orbit inherits its order from that of the group elements acting on the original starting vector (or ordered list of vectors.) When a given (set of) vector(s) arises more than once, the first instance defines the location of the vector(s) within the orbit.

mhi.mhi.make_momentum_spin_orbit(momenta, spin_dims, group, exchange_group=None)[source]#

Computes the momentum-spin orbit, i.e., the tensor product of the momentum orbit and the associated spinor orbit.

Parameters:
momenta(nmomenta, 3) array_like

The ordered momenta

spin_dimslist of ints

The dimensions of the spinor spaces.

group(|G|, 3, 3) array_like

The group matrices.

exchangearray_like or None

The exchange_group with namedtuple/WeightedPermutation elements

Returns:
spin_shellndarray

The flattened tensor product.

Notes

Consider a spinor transforming in an N-dimensional irrep. By definition, the different basis vectors for the irrep transform into each other under the action of the group. Thus, the spinorial part of the shell is just the tensor product of all the individual spinor spaces.

mhi.mhi.parity(permutation)[source]#

Computes the parity of permutation, assumed to be specified as a list of contiguous integers.

Parameters:
permutation: array_like

The permutation

Returns:
sign: +1 or -1

The parity.

Examples

>>> parity([2,3,4,5])
1
>>> parity([5,2,3,4])
-1
mhi.mhi.partition(arr)[source]#

Computes the partitions of the array according to the unique elements.

Parameters:
arrarray_like

The array to partition

Returns:
partitionsdict

The partitions, where the keys are the unique entries of the input array and the values are the associated indices in the input array.

Notes

The keys are sorted according to the first appearance in arr

Examples

>>> partition(['a','b','b','c']))
{'a': array([0]), 'b': array([1, 2]), 'c': array([3])}
mhi.mhi.recombine(labels, partition_keys, partition_idxs)[source]#

Recombines indices from partitioned sets. This function is equivalent to concatenating the indices for contiguous labels. For instance, [‘b1’, ‘b1’, ‘b2] are contiguous, but [‘b1’, ‘b2’, ‘b1’] are not.

Parameters:
labels(n, ) array_like

The labels, e.g., [‘b1’, ‘b2’, ‘b1’]

partition_keys(m, ) array_like

The unique elements of the labels, ordered by first appearance.

partition_idxs(m, ) tuple of array_like

The indices to recombine, in the same order as the partition_keys

Returns:
idxs(n, ) ndarray

The recombined indices

Examples

This example builds the tensor product of permutations on the set of labels [‘b1’,’b2’,’b1’]. Since the label ‘b1’ is not contiguous, the permutations associated with this label should include the identity (0)(2) and the swap (0,2).

>>> labels = ['b1','b2','b1']
>>> partitions = partition(labels)
>>> perms = {key: list(itertools.permutations(idxs)) for key, idxs in partitions.items()}
>>> keys = np.array(list(partitions.keys()))
>>> for idxs in itertools.product(*perms.values()):
>>>     print(mhi.recombine(labels, keys, idxs), "(correct)")
>>>     print(np.concatenate(idxs), "(wrong)")
[0 1 2] (correct)
[0 2 1] (wrong)
[2 1 0] (correct)
[2 0 1] (wrong)
mhi.mhi.make_exchange_projector(labels, tableau_map)[source]#

Computes a projector in the group algebra of the particle-exchange group.

Parameters:
labelslist

Particle labels, e.g., [‘a’,’b’,’c’,’b’,’c’,’b’]

tableau_mapdict

Young tableaux associated with each label, e.g., {‘a’: [[1]], ‘b’: [[1,2],[3]], ‘c’: [[1],[2]],}

Returns:
projectorlist of WeightedPermutation objects

The projection operator in the group algebra of the exchange group

mhi.mhi.make_exchange_projector_identical(labels)[source]#

Computes the exchange group associated with exchange of identical particles. The elements of this group are signed permuations.

Parameters:
labelsarray_like or None

The labels associated with the particles.

fermionsdict

Whether a given key (label) corresponds to a fermion. Values must be booleans: True for fermions and False for bosons.

Returns:
exchange_grouplist or None

The exchange group, with namedtuple/WeightedPermutation elements. The result is None when no labels are specified.

Notes

Generic fermions are labeled ‘f{n}’ for integers n, e.g., ‘f1’, ‘f2, … Generic bosons are labeld ‘b{n}’ for integers n, e.g., ‘b1’, ‘b2’, … Support also exists for certain named particles like ‘pi’, ‘K’, ‘neutron’, ‘proton’, ‘nucleon’.

mhi.mhi.make_internal_symmetry_projector(orbit, internal_symmetry)[source]#

Computes the projection matrix associated with an internal symmetry group.

Parameters:
momentum_spin_orbit(n, ) list of SpinShellTuple

Each element in the list defines the (momenta, spin) configuration for a given element of the orbit.

internal_grouplist of WeightedPermutation

The exchange group projector in the group algebra.

Returns:
proj(n, n) ndarray

projection matrix P, intended to be contracted against the extended momentum-spin representation matrices, giving projected matrices \(\hat{D}_{mm'}(R) = P_{mn} D_{nn'}(R) P_{nm'}\)

Notes

As a projection matrix, “proj” is idempotent: proj @ proj == proj.

The full internal symmetry operator can be thought of as a linear combination of permutations. The representation matrices can be computed either before or after taking the linear combination. This implementation takes the latter route, computing a representation matrix for each and taking a suitable linear combination.

mhi.mhi.multiply_perms(a, b)[source]#

Multiply the permutations a and b.

Group multiplication acts as composition when considering group elements as left-acting operators, i.e., \((a \cdot b)(v) = (a \circ b)(v) = a(b(v))\).

Parameters:
a, b(n,) ndarray
Returns:
ndarray

The product of permutations.

mhi.mhi.invert_perm(a)[source]#

Compute the inverse permutation such that \(a \circ a^{-1} = a^{-1} \circ a = (1)\)

Parameters:
a(n,) array_like
Returns:
ainv(n,) ndarray

The inverse permutation.

mhi.mhi.compose_permutation_algebra_elements(a, b)[source]#

Compute the product of two elements in the algebra of the permutation group.

Parameters:
a, blist of WeightedPermutation

The algebra elements to be multiplied.

Returns:
list

The product, as a list of WeightedPermutation objects.

mhi.mhi.algebra_elements_are_close(x1, x2)[source]#

Checks whether elements of the group algebra of Sn are numerically close.

Parameters:
x1, x2list of WeightedPermutation objects

The group algebra elements.

Returns:
are_closebool
mhi.mhi.is_valid_tableau(tableau)[source]#

Checks whether ragged list is a valid Young tableau.

Parameters:
tableaulist of lists

The candidate Young tableau.

Returns:
is_validbool
mhi.mhi.transpose_tableau(tableau)[source]#

Transposes a Young tableau.

Parameters:
tableaulist of lists

The tableau, given as a ragged list of integers.

Returns:
tableau_tlist of lists

The transposed tableau, as a ragged list of integers.

mhi.mhi.symmetrizer(row, *, signed, n)[source]#

Computes a symmetrizer over the specified indices.

The symmetrizer is an element of the algebra of the group Sn that yields the symmetrized or anti-symmetrized (depending on signed) space over the indices in row.

Parameters:
rowarray_like

The (1-indexed) indices over which to symmetrize.

signed0 or 1

Whether to symmetrize (0) or antisymmetrize (1)

nThe total number of indices
Returns:
list

The symmetrizer in the group algebra, represented as a list of WeightedPermutation objects.

Notes

Indices less than or equal to 1 are taken from the back of the list.

Examples

>>> mhi.symmetrizer([1,2], signed=1, n=3)
[WeightedPermutation(weight=1, perm=array([0, 1, 2])),
 WeightedPermutation(weight=-1, perm=array([1, 0, 2]))]
>>> mhi.symmetrizer([2,3], signed=0, n=3)
[WeightedPermutation(weight=1, perm=array([0, 1, 2])),
 WeightedPermutation(weight=1, perm=array([1, 0, 2]))]
>>> mhi.symmetrizer([0,-1], signed=0, n=3)
[WeightedPermutation(weight=1, perm=array([0, 1, 2])),
 WeightedPermutation(weight=1, perm=array([0, 2, 1]))]
mhi.mhi.make_young_projector(tableau, *, n)[source]#

Makes the projection operator, acting on n indices, associated with the given Young tableau.

Parameters:
tableaulist of lists

The Young tableau.

nint

The number of indices for which to construct the projector.

Returns:
projlist of WeightedPermutation objects

The projection operator in the group algebra.

Notes

Young tableaux are specified by indices 1,2,3,…. starting with 1. Indices in python arrays are zero indexed. This implementation uses definition of the Hermitian projection operators given by Ref. [1]. In particular, this function uses Eq (86) in Theorem 3 (KS Hermitian Young projectors). For large tableaux, it would likely advantageous to switch to the Measure Of Lexical Disorder (MOLD) definition of the projection operators given in Theorem 5.

References

[1]

J. Alcock-Zeilinger and H. Weigert “Compact Hermitian Young Projection Operators” J.Math.Phys. 58 (2017) 5, 051702 [arXiv:1610.10088].

mhi.mhi.make_momentum_orbit_rep_matrix(orbit, group_element)[source]#

Computes the momentum-orbit representation matrix associated with the action of a group element on an orbit of vectors.

Parameters:
orbitlist

The orbit of an ordered set of vectors under the group

group_elementarray_like

The group element which will act on all the elements of the orbit

Returns:
permutation(|O|, |O|) ndarray

The representation matrix, which happens to be a permutation matrix.

Notes

A group acts on a vector to generate an orbit O. When the group then acts on the orbit, the result is a permutation of the vectors in the orbit. The permutation can be represented as a square matrix of size |O|x|O|, where |O| is the size of the orbit. The full set of these matrices (with one for each group element) is itself a group representation.

mhi.mhi.make_momentum_orbit_rep(orbit, group)[source]#

Computes the representation of a group G acting on an orbit O.

Parameters:
orbitlist

The orbit of an ordered set of vectors under the group

group(|G|, 3, 3) array_like

The group matrices.

Returns:
representation(|G|, |O|, |O|) ndarray

The momentum-representation matrices \(D_{m,m'}(R)\) for all \(R \in G\).

mhi.mhi.make_momentum_spin_rep(Dmm, *Dspin)[source]#

Computes the combined momentum-spin representation matrices.

Parameters:
Dmm(|G|,|O|,|O|) ndarray

The momentum-representation matrices.

*Dspin(|G^D|, |\Gamma|, |\Gamma|) ndarray(s)

The spin irrep matrices.

Returns:
Dmm_spin(|G^D|, dim_total, dim_total) ndarray

The combined momomentum-spin representation matrices, where the total dimension is given by dim_total = |O|x|irrep1|x|irrep2|x...x|irrepN|. As should be expected, the proudct includes all the representations appearing in the list of “Dspin” matrices.

mhi.mhi.make_irrep_matrix(polarizations, group_element)[source]#

Computes the irrep matrix associated with a group element using the algebraic method in terms of polarization tensors.

Parameters:
polarizationslist

The polarization tensors specifying the basis functions for the irrep

group_element(3, 3) array_like

The group element.

Returns:
irrep_matrix(|\Gamma|, |\Gamma|) ndarray

The irrep matrix.

mhi.mhi.make_irrep(polarizations, group)[source]#

Computes the irrep matrices D_{mu, mu’}(g).

Parameters:
polarizationslist

The polarization tensors specifying the basis functions for the irrep.

group: ``(|G|, 3, 3)`` array_like

The group matrices.

Returns:
irrep(|G|, |\Gamma|, |\Gamma|) ndarray

The irrep matrices.

mhi.mhi.make_irrep_matrix_spinor(irrep_basis, group_element)[source]#

Computes the spinorial irrep matrix associated with a group element, given a basis for the irrep.

Parameters:
irrep_basislist

The spinor basis spanning the irrep

group_element(4, 4) ndarray

The “spinorial” group element.

Returns:
irrep_matrix(|\Gamma|, |\Gamma|) ndarray

The irrep matrix.

Notes

Each basis state should be a list of tuples “(coefficient, SpinorTuple)”

mhi.mhi.make_irrep_spinor(basis, group)[source]#

Computes the spinorial irrep matrices D_{mu, mu’}(g).

Parameters:
basislist

The spinor basis spanning the irrep

group(|G|, 4, 4) list or array_like

The group matrices.

Returns:
irrep(|G|, |\Gamma|, |\Gamma|) ndarray

The irrep matrices.

mhi.mhi.make_irrep_from_group(little_group)[source]#

Computes the irrep matrices associated with a particular little group.

Parameters:
little_group(|G|, 3, 3) ndarray

The matrices for the little group G.

Returns:
Dmumudict

The irrep matrices as a dict. The keys give the name of the irrep. The values contain the irrep matrices themselves, each with shape (|G|, |\Gamma|, |\Gamma|).

mhi.mhi.make_irrep_from_groupD(little_group)[source]#

Computes double-cover irrep matrices associated with a given little group, including both spinorial and bosonic irreps.

Parameters:
little_group(|G|, 3, 3) ndarray

The matrices for the little group G.

Returns:
Dmumu_doubledict

The irrep matrices as a dict. The keys give the name of the irrep. The values contain the irrep matrices themselves, each with shape (|G^D|, |\Gamma|, |\Gamma|).

mhi.mhi.orth(arr)[source]#

Computes an orthonormal basis for the row space.

This implementation constructs the basis for the row space using the Gram-Schmidt algorithm applied to the rows of the input array.

Parameters:
arrarray_like
Returns:
basisndarray

An orthonormal basis for the row space, with each row corresponding to a basis element.

mhi.mhi.project(vector, basis, direction)[source]#

Computes the paralell or perpendicular projection of a vector with respect to the space spanned by a basis.

Parameters:
vectorarray_like
basisarray_like

Vectors specifying the basis, with each row corresponding to a vector.

direction{‘parallel’, ‘perpendicular’}

Whether to compute the parallel or perpendicular projection.

Returns:
vector_newnp.array

The projected vector.

mhi.mhi.apply_schur(Dmm, Dmumu, verbose)[source]#

Computes the block diagonalization matrices using Schur’s algorithm, including transition operators to move between rows in a given irrep.

Parameters:
Dmm(|G|, |O|, |O|) array_like

Momentum-(spin)-representation matrices

Dmumudict

The irrep matrices. The keys give the name of the irrep. The values contain the group irrep matrices, with shape (|G|, |\Gamma|, |\Gamma|).

Returns:
u_matrixdict

The block diagonalization matrix. The keys are tuples of the form (irrep_name, degeneracy_number). The values are arrays containing the matrices, each of shape (|\Gamma|, |O|).

transition_operators: dict

Column of transition operators T_{\mu,0}. The keys are irrep names. The values are ndarrays of shape (|\Gamma|, |O|, |O|).

mhi.mhi.rephase(arr, irrep)[source]#

Applies a phase convention to block diagonalization matrices

Parameters:
arr(|\Gamma|, |O|) array_like

Table specifiying the matrix

Returns:
arr(|\Gamma|, |O|) ndarray

The table with appropriate phases applied.

Notes

The phase convention is as follows:
  • For a generic irrep, the phase is chosen such that the first nonzero entry of the first row (\(\mu=1\)) is real and positive.

  • For the irreps \(T_2^+\) and \(T_2^-\) only, the phase is chosen such that the second row (\(\mu=2\)) is purely imaginary with a negative imaginary part. This choice matches the basis-vector conventions of Basak et al., where a particular combination of spheric harmonics \((Y_2^2 - Y_2^{-2})\) is used as the \(\mu=2\) basis vector for T2.

References

[1]

S. Basak et al., “Clebsch-Gordan construction of lattice interpolating fields for excited baryons”, Phys. Rev. D 72, 074501 (2005), [arXiv:hep-lat/0508018].

mhi.mhi.load_particle_info(fname=None)[source]#

Loads tabulated information for particle names, spins (boson vs fermion), and spin irreps.

Parameters:
fnamestr or None

The path to the input yaml file with the tabulated data

Returns:
tabledict

The particle information in the form {<name> : (<fermion?>, <irrep>}

mhi.mhi.make_pseudoscalar_irrep(little)[source]#

Instantiates the double-cover irrep matrices associated with a pseudoscalar particle, assumed to transform under the \(A_1^-\) irrep. Usually this function is used when constructing the Dspin matrices.

Parameters:
little(|G|, 3, 3) array_like

The little-group matrices associated with some momenta

Returns:
pseudoscalar(|G^D|, 1, 1) ndarray

The A1m irrep matrices restricted to the double cover of the little group.

mhi.mhi.make_spin_half_irrep(little_double)[source]#

Instantiates the double-cover irrep matrices associated with a spin-half particle, assumed to transform under the \(G_1^+\) irrep. Usually this function is used when constructing the Dspin matrices.

Parameters:
little_double(|G^D|, 4, 4) array_like

The spinorial little-group matrices associated with some momenta

Returns:
pseudoscalar(|G^D|, 1, 1) ndarray

The G_1^+ irrep matrices restricted to the double cover of the little group.

mhi.mhi.make_Dspin(spin_irreps, little, little_double)[source]#

Builds a list of spin irrep matrices for the specified particles.

Parameters:
spin_irrepslist

The names of the particles’ spin irreps as strings, e.g., [‘A1m’, ‘G1p’]

little(|G|, 3, 3) array_like

The little-group matrices associated with some momenta

little_double(|G^D|, 4, 4) array_like

The spinorial little-group matrices associated with some momenta

Returns:
Dspinlist

The “spin” irrep matrices for each particle

mhi.mhi.identify_spin_dim(irrep)[source]#

Identifies the dimension of the specified spin irrep.

Parameters:
strirrep

The name of the spin irrep, e.g., ‘A1m’ or ‘G1p’.

Returns:
dimint

The dimension of the irrep

class mhi.mhi.IrrepDecomposition(decomp, orbit, Dmm, Dmumu, little_name, stab_name, transition_operators)[source]#

Container for results of computing the block-diagonalization matrices which project linear combinations of plane-wave states onto irreps of the cubic group.

Parameters:
decompdict

The irrep decomposition and change-of-basis matrices. The keys are tuples (irrep_name, degeneracy_idx). The values are the block-diagonalization matrices, given as arrays of shape (|Gamma|, |O|).

orbit(|O|,) list of SpinShellTuple

The “extended” spin-momentum orbit, where each element is a SpinShellTuple specifying momentum and spin indices.

Dmm(|G|, |O|, |O|), ndarray

The (reducible) representation matrices associated with the orbit.

Dmumu(|G|, |Gamma|, |Gamma|), ndarray

The irrep matrices

little_namestr

The name of the little group leaving the total momentum invariant

stab_namestr

The name of the stabilizer group leaving the ordered set of momenta invariant.

transition_operators: dict

Column of transition operators \(T_{\mu,0}\). The keys are irrep names. The values are ndarrays of shape (|Gamma|, |O|, |O|).

Notes

Projection onto cubic-group irreps requires two pieces.
  • A basis of momentum plane-wave correlation functions, presumably computed using lattice QCD.

  • The block-diagonalization matrices computed using this module.

To carry out this projection, the basis of correlation functions must have the same order as the orbit used to compute the block-diagonalization matrices. The required ordering can be seen by examining the “orbit.”

format(latex=False)[source]#

Formats the irrep decomposition as text string.

Parameters:
latexbool

Whether or not return a latex-formatted sting

Returns:
outputstr

The irrep decomposition, e.g., "A1p + A2p + 2*Ep" or "$A_1^+ \oplus A_2^+ \oplus 2E^+$".

little_name(latex=False)[source]#

Returns the name of the little group.

Parameters:
latexbool

Whether or not return a latex-formatted sting

Returns:
little_namestr

The name of the little group

stab_name(latex=False)[source]#

Returns the name of the stabilizer group.

Parameters:
latexbool

Whether or not return a latex-formatted sting

Returns:
stab_namestr

The name of the stabilizer group

mhi.mhi.test_stabilizer(verbose=False)[source]#

Test identification of the stabilizer group by comparing to known results.

mhi.mhi.anticommutator(arr1, arr2)[source]#

Computes the anticommutator of two matrices.

Parameters:
arr1array_like
arr2array_like
Returns:
arrndarray

The anticommutator {arr1, arr2}

mhi.mhi.commutator(arr1, arr2)[source]#

Computes the commutator of two matrices.

Parameters:
arr1array_like
arr2array_like
Returns:
arrndarray

The commutator {arr1, arr2}

mhi.mhi.test_clifford(gamma, eta, verbose=False)[source]#

Tests the Clifford-algebra condition, \(\{\gamma_\mu, \gamma_\nu\} = 2*\eta_{\mu\nu}\)

Parameters:
gamma(4, ) list or (4, 4, 4) array_like

The gamma matrices gamma[i], i=0,1,2,3.

eta(4,4) ndarray

The metric.

verbosebool

Whether to print additional information about successful tests.

Returns:
None
mhi.mhi.test_gamma5(gamma, gamma5, verbose)[source]#

Tests anticommutation of gamma5, \(\{\gamma_\mu, \gamma_5\} = 0\)

Parameters:
gamma(4, ) list or (4, 4, 4) array_like

The gamma matrices gamma[i], i=0,1,2,3.

gamma5(4,4)

The matrix gamma5.

verbosebool

Whether to print additional information about successful tests.

Returns:
None
mhi.mhi.test_row_orthogonality(u_matrix, verbose=False)[source]#

Tests the row orthogonality of tables of block-diagonalization matrices, for a given irrep.

Parameters:
u_matrixdict

The block diagonalization matrix. The keys are tuples of the form (irrep_name, degeneracy_number). The values are arrays containing the matrices, each of shape (|Gamma|, |O|).

verbosebool

Whether to print additional information about successful tests.

Returns:
None
mhi.mhi.count_degeneracy(u_matrix)[source]#

Counts the degeneracy of irreps within the tables of block-diagnoalization matrices.

Parameters:
u_matrixdict

The block diagonalization matrix. The keys are tuples of the form (irrep_name, degeneracy_number). The values are arrays containing the matrices, each of shape (|Gamma|, |O|).

Returns:
degeneracydict

The degeneracies, with keys corresponding to the irrep names and values to the counts.

mhi.mhi.test_degenerate_orthogonality(u_matrix, verbose=False)[source]#

Tests that tables corresponding to degenerate irreps are orthogonal.

Parameters:
u_matrixdict

The block diagonalization matrix. The keys are tuples of the form (irrep_name, degeneracy_number). The values are arrays containing the matrices, each of shape (|Gamma|, |O|).

verbosebool

Whether to print additional information about successful tests.

Returns:
None
mhi.mhi.test_block_diagonalization(Dmm, Dmumu, u_matrix, verbose=False)[source]#

Tests the block diagonalization property of the change-of-basis matrices, \(D_{\mu\nu} = U^\dagger_{\mu m} D_{mm'} U_{m' \nu}\) in terms of the given momentum-shell representation Dmm, the block-diagonalization matrix U, and the block-diagonal irrep matrix Dmumu.

Parameters:
Dmm(|G|,|O|,|O|) ndarray

The momentum-representation matrices.

Dmumudict

The irrep matrices as a dict. The keys give the name of the irrep. The values contain the irrep matrices themselves, each with shape (|G|, |Gamma|, |Gamma|).

u_matrixdict

The block diagonalization matrix. The keys are tuples of the form (irrep_name, degeneracy_number). The values are arrays containing the matrix, each of shape (|Gamma|, |O|).

verbosebool

Whether to print additional information about successful tests.

Returns:
None
mhi.mhi.write_hdf5(h5fname, result_dict)[source]#

Writes block-diagonalization matrices to HDF5.

Parameters:
h5fnamestr

The name of the output file

result_dictdict

The block-diagonalization / change of basis matrices, where each key is a tuple the form (irrep_name, degeneracy) and each key is an ndarray with the associated matrix.

Returns:
None
mhi.mhi.read_hdf5(h5fname)[source]#

Reads saved block-diagonalization matrices from HDF5.

Parameters:
h5fnamestr

The name of HDF5 file to read

Returns:
result_dictdict

The block-diagonalization / change of basis matrices, where each key is a tuple the form (irrep_name, degeneracy) and each key is an ndarray with the associated matrix.

Basis functions module#

class mhi.basis_functions.SpinorTuple(j, jz, parity)[source]#

A spinor basis state |j, jz, parity>.

Attributes:
jfloat

The half-integer total spin.

jzfloat

The half-integer z component of spin.

parity+1 or -1

The parity.

mhi.basis_functions.basis_Oh_A1p = [0.577350269189626*x**2 + 0.577350269189626*y**2 + 0.577350269189626*z**2]#

Basis functions for irrep \(A_{1}^{+}\) of \(O_h\).

mhi.basis_functions.basis_Oh_A1m = [0.408248290463863*x*y*z*(x**4*(y**2 - z**2) + y**4*(-x**2 + z**2) + z**4*(x**2 - y**2))]#

Basis functions for irrep \(A_{1}^{-}\) of \(O_h\).

mhi.basis_functions.basis_Oh_A2p = [0.408248290463863*x**4*(y**2 - z**2) + 0.408248290463863*y**4*(-x**2 + z**2) + 0.408248290463863*z**4*(x**2 - y**2)]#

Basis functions for irrep \(A_{2}^{+}\) of \(O_h\).

mhi.basis_functions.basis_Oh_A2m = [x*y*z]#

Basis functions for irrep \(A_{2}^{-}\) of \(O_h\).

mhi.basis_functions.basis_Oh_Ep = [-0.408248290463863*x**2 - 0.408248290463863*y**2 + 0.816496580927726*z**2, 0.707106781186547*x**2 - 0.707106781186547*y**2]#

Basis functions for irrep \(E^{+}\) of \(O_h\).

mhi.basis_functions.basis_Oh_Em = [x*y*z*(x**2 - y**2), -x*y*z*(-x**2 - y**2 + 2*z**2)]#

Basis functions for irrep \(E^{-}\) of \(O_h\).

mhi.basis_functions.basis_Oh_T1m = [1.4142135623731*z, -x + 1.0*I*y, x + 1.0*I*y]#

Basis functions for irrep \(T_1^-\) of \(O_h\).

mhi.basis_functions.basis_Oh_T1p = [1.4142135623731*x*y*(x**2 - y**2), 1.0*I*x*z*(-x**2 + z**2) - y*z*(y**2 - z**2), 1.0*I*x*z*(-x**2 + z**2) + y*z*(y**2 - z**2)]#

Basis functions for irrep \(T_1^+\) of \(O_h\).

mhi.basis_functions.basis_Oh_T2m = [1.0*I*x*(y**2 - z**2) - y*(-x**2 + z**2), -1.4142135623731*I*z*(x**2 - y**2), 1.0*I*x*(y**2 - z**2) + y*(-x**2 + z**2)]#

Basis functions for irrep \(T_2^-\) of \(O_h\).

mhi.basis_functions.basis_Oh_T2p = [-x*z + 1.0*I*y*z, -1.4142135623731*I*x*y, x*z + 1.0*I*y*z]#

Basis functions for irrep \(T_2^+\) of \(O_h\).

mhi.basis_functions.basis_C4v_A1 = [0.577350269189626*x**2 + 0.577350269189626*y**2 + 0.577350269189626*z**2]#

Basis functions for irrep \(A_{1}\) of \(C_{4v}\) = \(A_{1}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C4v_A2 = [0.408248290463863*x*y*z*(x**4*(y**2 - z**2) + y**4*(-x**2 + z**2) + z**4*(x**2 - y**2))]#

Basis functions for irrep \(A_{2}\) of \(C_{4v}\) = \(A_{1}^{-}\) of \(O_h\).

mhi.basis_functions.basis_C4v_B1 = [0.707106781186547*x**2 - 0.707106781186547*y**2]#

Basis functions for irrep \(B_{1}\) of \(C_{4v}\) = \(E_{1}^{+} (\mu=2)\) of \(O_h\).

mhi.basis_functions.basis_C4v_B2 = [-1.4142135623731*I*x*y]#

Basis functions for irrep \(B_{2}\) of \(C_{4v}\) = \(T_2^+ (\mu=2)\) of \(O_h\).

mhi.basis_functions.basis_C4v_E = [-x*z, -y*z]#

Basis functions for irrep \(E\) of \(C_{4v}\).

mhi.basis_functions.basis_C3v_A1 = [0.577350269189626*x**2 + 0.577350269189626*y**2 + 0.577350269189626*z**2]#

Basis functions for irrep \(A_{1}\) of \(C_{3v}\) = \(A_{1}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C3v_A2 = [0.408248290463863*x**4*(y**2 - z**2) + 0.408248290463863*y**4*(-x**2 + z**2) + 0.408248290463863*z**4*(x**2 - y**2)]#

Basis functions for irrep \(A_{2}\) of \(C_{3v}\) = \(A_{2}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C3v_E = [-0.408248290463863*x**2 - 0.408248290463863*y**2 + 0.816496580927726*z**2, 0.707106781186547*x**2 - 0.707106781186547*y**2]#

Basis functions for irrep \(E\) of \(C_{3v}\) = \(E_{}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C2v_A1 = [0.577350269189626*x**2 + 0.577350269189626*y**2 + 0.577350269189626*z**2]#

Basis functions for irrep \(A_{1}\) of \(C_{2v}\) = \(A_{1}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C2v_A2 = [0.408248290463863*x*y*z*(x**4*(y**2 - z**2) + y**4*(-x**2 + z**2) + z**4*(x**2 - y**2))]#

Basis functions for irrep \(A_{2}\) of \(C_{2v}\) = \(A_{1}^{-}\) of \(O_h\).

mhi.basis_functions.basis_C2v_B1 = [x*y*z]#

Basis functions for irrep \(B_{1}\) of \(C_{2v}\) = \(A_{2}^{-}\) of \(O_h\).

mhi.basis_functions.basis_C2v_B2 = [0.408248290463863*x**4*(y**2 - z**2) + 0.408248290463863*y**4*(-x**2 + z**2) + 0.408248290463863*z**4*(x**2 - y**2)]#

Basis functions for irrep \(B_{2}\) of \(C_{2v}\) = \(A_{2}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C2R_A = [0.577350269189626*x**2 + 0.577350269189626*y**2 + 0.577350269189626*z**2]#

Basis functions for irrep \(A\) of \(C_2^R\) = \(A_{1}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C2R_B = [0.408248290463863*x*y*z*(x**4*(y**2 - z**2) + y**4*(-x**2 + z**2) + z**4*(x**2 - y**2))]#

Basis functions for irrep \(B\) of \(C_2^R\) = \(A_{1}^{-}\) of \(O_h\).

mhi.basis_functions.basis_C2P_A = [0.577350269189626*x**2 + 0.577350269189626*y**2 + 0.577350269189626*z**2]#

Basis functions for irrep \(A\) of \(C_2^P\) = \(A_{1}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C2P_B = [0.408248290463863*x**4*(y**2 - z**2) + 0.408248290463863*y**4*(-x**2 + z**2) + 0.408248290463863*z**4*(x**2 - y**2)]#

Basis functions for irrep \(B\) of \(C_2^P\) = \(A_{2}^{+}\) of \(O_h\).

mhi.basis_functions.basis_C1_A = [0.577350269189626*x**2 + 0.577350269189626*y**2 + 0.577350269189626*z**2]#

Basis functions for irrep \(A\) of \(C_1\) = \(A_{1}^{+}\) of \(O_h\).

mhi.basis_functions.basis_OhD_G1p = [[(1, (0.5, 0.5, 1))], [(1, (0.5, -0.5, 1))]]#

Basis spinors for irrep \(G_{1}^{+}\) of \(O_h^D\).

mhi.basis_functions.basis_OhD_G1m = [[(1, (0.5, 0.5, -1))], [(1, (0.5, -0.5, -1))]]#

Basis spinors for irrep \(G_{1}^{-}\) of \(O_h^D\).

mhi.basis_functions.basis_OhD_G2p = [[(0.408248290463863, (2.5, 2.5, 1)), (-0.9128709291752769, (2.5, -1.5, 1))], [(0.408248290463863, (2.5, -2.5, 1)), (-0.9128709291752769, (2.5, 1.5, 1))]]#

Basis spinors for irrep \(G_{2}^{+}\) of \(O_h^D\).

mhi.basis_functions.basis_OhD_G2m = [[(0.408248290463863, (2.5, 2.5, -1)), (-0.9128709291752769, (2.5, -1.5, -1))], [(0.408248290463863, (2.5, -2.5, -1)), (-0.9128709291752769, (2.5, 1.5, -1))]]#

Basis spinors for irrep \(G_{2}^{-}\) of \(O_h^D\).

mhi.basis_functions.basis_OhD_Hp = [[(1, (1.5, 1.5, 1))], [(1, (1.5, 0.5, 1))], [(1, (1.5, -0.5, 1))], [(1, (1.5, -1.5, 1))]]#

Basis spinors for irrep \(H^{+}\) of \(O_h^D\).

mhi.basis_functions.basis_OhD_Hm = [[(1, (1.5, 1.5, -1))], [(1, (1.5, 0.5, -1))], [(1, (1.5, -0.5, -1))], [(1, (1.5, -1.5, -1))]]#

Basis spinors for irrep \(H^{-}\) of \(O_h^D\).

mhi.basis_functions.basis_Dic4_G1 = [[(1, (0.5, 0.5, 1))], [(1, (0.5, -0.5, 1))]]#

Basis spinors for irrep \(G_{1}\) of \(\mathrm{Dic}_4\).

mhi.basis_functions.basis_Dic4_G2 = [[(0.408248290463863, (2.5, 2.5, 1)), (-0.9128709291752769, (2.5, -1.5, 1))], [(0.408248290463863, (2.5, -2.5, 1)), (-0.9128709291752769, (2.5, 1.5, 1))]]#

Basis spinors for irrep \(G_{2}\) of \(\mathrm{Dic}_4\).

mhi.basis_functions.basis_Dic3_G = [[(1, (0.5, 0.5, 1))], [(1, (0.5, -0.5, 1))]]#

Basis spinors for irrep \(G\) of \(\mathrm{Dic}_3\).

mhi.basis_functions.basis_Dic3_F1 = [[(0.5, (1.5, 1.5, 1)), ((0.08455098936288137+0.49279927982674443j), (1.5, 0.5, 1)), ((0.4082482904638631+0.2886751345948129j), (1.5, -0.5, 1)), ((-0.35355339059327373-0.35355339059327373j), (1.5, -1.5, 1))]]#

Basis spinors for irrep \(F_{1}\) of \(\mathrm{Dic}_3\).

mhi.basis_functions.basis_Dic3_F2 = [[(0.5, (1.5, 1.5, 1)), ((0.4927992798267444+0.08455098936288137j), (1.5, 0.5, 1)), ((-0.4082482904638631+0.2886751345948129j), (1.5, -0.5, 1)), ((0.35355339059327373+0.35355339059327373j), (1.5, -1.5, 1))]]#

Basis spinors for irrep \(F_{2}\) of \(\mathrm{Dic}_3\).

mhi.basis_functions.basis_Dic2_G = [[(1, (0.5, 0.5, 1))], [(1, (0.5, -0.5, 1))]]#

Basis spinors for irrep \(G\) of \(\mathrm{Dic}_2\).

mhi.basis_functions.basis_C4R_F1 = [[(1, (0.5, 0.5, 1))]]#

Basis spinors for irrep \(F_{1}\) of \(C_4^R\).

mhi.basis_functions.basis_C4R_F2 = [[(1, (0.5, -0.5, 1))]]#

Basis spinors for irrep \(F_{2}\) of \(C_4^R\).

mhi.basis_functions.basis_C4P_F1 = [[(0.7071067811865475, (0.5, 0.5, 1)), ((0.5-0.5j), (0.5, -0.5, 1))]]#

Basis spinors for irrep \(F_{1}\) of \(C_4^P\).

mhi.basis_functions.basis_C4P_F2 = [[(0.7071067811865475, (0.5, 0.5, 1)), ((-0.5+0.5j), (0.5, -0.5, 1))]]#

Basis spinors for irrep \(F_{2}\) of \(C_4^P\).

mhi.basis_functions.basis_C1D_F = [[(1, (0.5, 0.5, 1))]]#

Basis spinors for irrep \(F\) of \(C_1^D\).

mhi.basis_functions.basis_nucleon = [[(1, (0.5, 0.5, 1))], [(1, (0.5, -0.5, 1))]]#

Basis spinors for the nucleon.