Isotropy subgroup#

class spgrep_modulation.isotropy.IsotropyEnumerator(little_rotations: numpy.ndarray[Any, numpy.dtype[numpy.int32]], little_translations: numpy.ndarray[Any, numpy.dtype[numpy.float64]], qpoint: numpy.ndarray[Any, numpy.dtype[numpy.float64]], small_rep: numpy.ndarray[Any, numpy.dtype[numpy.complex128]], max_denominator: int = 100, atol: float = 1e-06)[source]#

Enumerate isotropy subgroups of given little group and small representation.

Initialize class.

Parameters
  • little_rotations (array, (order, 3, 3)) –

  • little_translations (array, (order, 3)) –

  • qpoint (array, (3, )) –

  • small_rep ((order, dim, dim)) – Small representation of space group at qpoint

  • max_denominator (int) –

  • atol (float) –

property little_rotations: numpy.ndarray[Any, numpy.dtype[numpy.int32]]#

Return rotations of little group at qpoint.

property little_translations: numpy.ndarray[Any, numpy.dtype[numpy.float64]]#

Return translations of little group at qpoint.

property qpoint: numpy.ndarray[Any, numpy.dtype[numpy.float64]]#

Return qpoint.

property small_rep: numpy.ndarray[Any, numpy.dtype[numpy.complex128]]#

Return small representation of space group.

property maximal_isotropy_subgroups: list[list[int]]#

Return list of indices for isotropy subgroups.

Let subgroup = self.maximal_isotropy_subgroups[i]. (self.little_rotations[subgroup], self.little_translations[subgroup]) gives a coset of the i-th isotropy subgroup.

property order_parameter_directions: list[numpy.ndarray[Any, numpy.dtype[numpy.complex128]]]#

Return order-parameter directions of isotropy subgroups.

Let opd = self.order_parameter_directions[i], which is an array with (num_direction, dim). num_direction is the number of free parameters for the i-th isotropy subgroup. dim is dimension of the given small representation.

spgrep_modulation.isotropy.get_translational_subgroup(qpoint: list[float] | NDArrayFloat, max_denominator: int = 100)[source]#

Return transformation matrix of the following sublattice.

Let t be a lattice point of the returned sublattice. Then, np.dot(t, qpoint) is integer.

Returns

transformationtransformation @ qpoint is integer vector.

Return type

array, (3, 3)

spgrep_modulation.isotropy.enumerate_point_subgroup(table: numpy.ndarray[Any, numpy.dtype[numpy.int32]], preserve_sublattice: list[bool], return_conjugacy_class: bool = True) list[list[int]][source]#

Enumerate conjugacy subgroups of point group.

spgrep_modulation.isotropy.search_compliment(rotations, translations, indices: list[int], transformation: numpy.ndarray[Any, numpy.dtype[numpy.int32]], table: numpy.ndarray[Any, numpy.dtype[numpy.int32]], atol: float = 1e-06) bool[source]#

Return true if some adjusted translations enable given point group and translational subgroup correspond to a space group.