spacegroup
- class irrep.spacegroup.SpaceGroup(cell, spinor=True, refUC=None, shiftUC=None, search_cell=False, trans_thresh=1e-05, alat=None, from_sym_file=None, v=0)[source]
Determine the space-group and save info in attributes. Contains methods to describe and print info about the space-group.
- Parameters:
cell (tuple, default=None) – cell[0] is a 3x3 array where cartesian coordinates of basis vectors a, b and c are given in rows. cell[1] is an array where each row contains the direct coordinates of an ion’s position. cell[2] is an array where each element is a number identifying the atomic species of an ion. See cell parameter of function get_symmetry in Spglib.
spinor (bool, default=True) – True if wave-functions are spinors (SOC), False if they are scalars.
refUC (array, default=None) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
shiftUC (array, default=None) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
search_cell (bool, default=False) – Whether the transformation to the conventional cell should be computed. It is True if kpnames was specified in CLI.
trans_thresh (float, default=1e-5) – Threshold used to compare translational parts of symmetries.
alat (float, default=None) – Lattice parameter in angstroms (quantum espresso convention).
from_sym_file (str, default=None) – If provided, the symmetry operations are read from this file. (format of pw2wannier90 prefix.sym file)
v (int, default=0) – Verbosity level. Default set to minimalistic printing
- spinor
True if wave-functions are spinors (SOC), False if they are scalars.
- Type:
bool
- symmetries
Each element is an instance of class SymmetryOperation corresponding to a symmetry in the point group of the space-group.
- Type:
list
- symmetries_tables
Attribute symmetries of class IrrepTable. Each component is an instance of class SymopTable corresponding to a symmetry operation in the “point-group” of the space-group.
- Type:
list
- name
Symbol of the space-group in Hermann-Mauguin notation.
- Type:
str
- number
Number of the space-group.
- Type:
int
- Lattice
Each row contains cartesian coordinates of a basis vector forming the unit-cell in real space.
- Type:
array, shape=(3,3)
- positions
Direct coordinate of sites in the DFT cell setting.
- Type:
array
- typat
Indices to identify the element in each atom. Atoms of the same element share the same index.
- Type:
list
- order
Number of symmetries in the space group (in the coset decomposition w.r.t. the translation subgroup).
- Type:
int
- refUC
3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
- Type:
array, default=None
- shiftUC
Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
- Type:
array, default=None
- alat
Lattice parameter in angstroms (quantum espresso convention).
- Type:
float
- from_sym_file
if provided, the symmetry operations are read from this file. (format of pw2wannier90 prefix.sym file)
- Type:
str, default=None
Notes
The symmetry operations to which elements in the attribute symmetries correspond are the operations belonging to the point group of the space-group \(G\), i.e. the coset representations \(g_i\) taking part in the coset decomposition of \(G\) w.r.t the translation subgroup \(T\): ..math:: G=T+ g_1 T + g_2 T +…+ g_N T
- determine_basis_transf(refUC_cli, shiftUC_cli, refUC_lib, shiftUC_lib, search_cell, trans_thresh, v=0)[source]
Determine basis transformation to conventional cell. Priority is given to the transformation set by the user in CLI.
- Parameters:
refUC_cli (array) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting. Set in CLI.
shiftUC_cli (array) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting. Set in CLI.
refUC_lib (array) – Obtained via spglib.
shiftUC_lib (array) – Obtained via spglib. It may not be the shift taking the origin to the position adopted in the tables (BCS). For example, origin choice 1 of ITA is adopted in spglib for centrosymmetric groups, while origin choice 2 in BCS.
search_cell (bool, default=False) – Whether the transformation to the conventional cell should be computed. It is True if kpnames was specified in CLI.
trans_thresh (float, default=1e-5) – Threshold to compare translational parts of symmetries.
v (int, default=0) – Verbosity level. Default set to minimalistic printing
- Returns:
array – Transformation of vectors defining the unit cell to the standard setting.
array – Shift taking the origin of the unit cell used in the DFT calculation to that of the convenctional setting used in BCS.
- Raises:
RuntimeError – Could not find a pait (refUC,shiftUC) matching symmetries obtained from spglib to those in tables (mod. lattice translations of the primitive cell).
- get_irreps_from_table(kpname, K, v=0)[source]
Read irreps of the little-group of a maximal k-point.
- Parameters:
kpname (str) – Label of the maximal k-point.
K (array, shape=(3,)) – Direct coordinates of the k-point.
v (int, default=0) – Verbosity level. Default set to minimalistic printing
- Returns:
tab – Each key is the label of an irrep, each value another dict. Keys of every secondary dict are indices of symmetries (starting from 1 and following order of operations in tables of BCS) and values are traces of symmetries.
- Return type:
dict
- Raises:
RuntimeError – Translational or rotational parts read from tables and found by spglib do not match for a symmetry operation.
RuntimeError – A symmetry from the tables matches with many symmetries found by spglib.
RuntimeError – The label of a k-point given in the CLI matches with the label of a k-point read from tables, but direct coordinates of these k-vectors do not match.
RuntimeError – There is not any k-point in the tables whose label matches that given in parameter kpname.
- json(symmetries=None)[source]
Prepare dictionary with info of space group to save in JSON
- Returns:
d – Dictionary with info about space group
- Return type:
dict
- match_symmetries(refUC=None, shiftUC=None, signs=False, trans_thresh=1e-05)[source]
Matches symmetry operations of two lists. Translational parts are matched mod. lattice translations (important for centered structures).
- Parameters:
refUC (array, default=None) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
shiftUC (array, default=None) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
signs (bool, default=False) – If True, match also rotations of spinors corresponding to each symmetry.
trans_thresh (float, default=1e-5) – Threshold used to compare translational parts of symmetries.
- Returns:
list – The \(i^{th}\) element corresponds to the \(i^{th}\) symmetry found by spglib and it is the position that the same symmetry has in the tables.
list – The \(i^{th}\) element corresponds to the \(i^{th}\) symmetry found by spglib and it is the phase difference w.r.t. the same symmetry in the tables, which may arise if their translational parts differ by a lattice vector.
array – The \(i^{th}\) element corresponds to the \(i^{th}\) symmetry found by spglib and it is the sign needed to make the matrix for spinor rotation identical to that in tables.
Note
Arguments symmetries1 and symmetries2 always take values of attributes self.symmetry and self.symmetries_tables, so they can be removed, but this way keeps the function more generic.
- show(symmetries=None)[source]
Print description of space-group and symmetry operations.
- Parameters:
symmetries (int, default=None) – Index of symmetry operations whose description will be printed. Run IrRep with flag onlysym to check the index corresponding to each symmetry operation.
- str()[source]
Create a string to describe of space-group and its symmetry operations.
- Returns:
Description to print.
- Return type:
str
- vecs_centering()[source]
Check the space group and generate vectors of centering.
- Returns:
Each row is a lattice vector describing the centering.
- Return type:
array
- vecs_inv_centers()[source]
Get the positions of all inversion centers in the unit cell.
- Returns:
Each element is a vector pointing to a position center in the convenctional unit cell.
- Return type:
array
Notes
If the space group is primitive, there are 8 inversion centers, but there are more if it is a group with a centering.
- write_sym_file(filename, alat=None)[source]
Write symmetry operations to a file.
- Parameters:
filename (str) – Name of the file.
alat (float, default=None) – Lattice parameter in angstroms. If not specified, the lattice parameter is not written to the file.
- write_trace()[source]
Construct description of matrices of symmetry operations of the space-group in the format: {{R|t}}-> R11,R12,…,R23,R33,t1,t2,t3 and, when SOC was included, the elements of the matrix describing the transformation of the spinor in the format: Re(S11),Im(S11),Re(S12),…,Re(S22),Im(S22).
- Returns:
String describing matrices of symmetry operations.
- Return type:
str
- class irrep.spacegroup.SymmetryOperation(rot, trans, Lattice, ind=-1, spinor=True)[source]
Contains information to describe a symmetry operation and methods to get info about the symmetry, transform it to a reference choice of unit cell and print a description of it.
- Parameters:
rotation (array, shape=(3,3)) – Matrix describing the tranformation of basis vectors of the unit cell under the symmetry operation.
translation (array, shape=(3,)) – Translational part of the symmetry operation, in terms of the basis vectors of the unit cell.
Lattice (array, shape=(3,3)) – Each row contains cartesian coordinates of a basis vector forming the unit-cell in real space.
ind (int, default=-1) – Index of the symmetry operation.
spinor (bool, default=true) – True if wave-functions are spinors, False if they are scalars.
- ind
Index of the symmetry operation.
- Type:
int
- rotation
Matrix describing the tranformation of basis vectors of the unit cell under the symmetry operation.
- Type:
array, shape=(3,3)
- translation
Translational part of the symmetry operation, in terms of the basis vectors of the unit cell.
- Type:
array, shape=(3,)
- Lattice
Each row contains cartesian coordinates of a basis vector forming the unit-cell in real space.
- Type:
array, shape=(3,3)
- axis
Rotation axis of the symmetry.
- Type:
array, shape=(3,)
- angle
Rotation angle of the symmmetry, in radians.
- Type:
float
- inversion
False if the symmetry preserves handedness (identity, rotation, translation or screw rotation), True otherwise (inversion, reflection roto-inversion or glide reflection).
- Type:
bool
- angle_str
String describing the rotation angle in radians.
- Type:
str
- spinor
True if wave-functions are spinors, False if they are scalars.
- Type:
bool
- spinor_rotation
Matrix describing how spinors transform under the symmetry.
- Type:
array, shape=(2,2)
- sign
Factor needed to match the matrix for the rotation of spinors to that in tables.
- Type:
float
- get_angle_str()[source]
Give str of rotation angle.
- Returns:
Rotation angle in radians.
- Return type:
str
- Raises:
RuntimeError – Angle does not belong to 1, 2, 3, 4 or 6-fold rotation.
- json_dict(refUC=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), shiftUC=array([0., 0., 0.]))[source]
Prepare dictionary with info of symmetry to save in JSON
- Returns:
d – Dictionary with info about symmetry
- Return type:
dict
- rotation_refUC(refUC)[source]
Calculate the matrix of the symmetry in the reference cell choice.
- Parameters:
refUC (array) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
- Returns:
R1 – Matrix for the transformation of basis vectors forming the reference unit cell.
- Return type:
array, shape=(3,3)
- Raises:
RuntimeError – If the matrix contains non-integer elements after the transformation.
- show(refUC=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), shiftUC=array([0., 0., 0.]))[source]
Print description of symmetry operation.
- Parameters:
refUC (array, default=np.eye(3)) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
shiftUC (array, default=np.zeros(3)) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
- str(refUC=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), shiftUC=array([0., 0., 0.]))[source]
Construct description of symmetry operation.
- Parameters:
refUC (array, default=np.eye(3)) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
shiftUC (array, default=np.zeros(3)) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
- Returns:
Description to print.
- Return type:
str
- str2(refUC=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), shiftUC=array([0., 0., 0.]))[source]
Print matrix of a symmetry operation in the format: {{R|t}}-> R11,R12,…,R23,R33,t1,t2,t3 and, when SOC was included, the elements of the matrix describing the transformation of the spinor in the format: Re(S11),Im(S11),Re(S12),…,Re(S22),Im(S22).
- Parameters:
refUC (array, default=np.eye(3)) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
shiftUC (array, default=np.zeros(3)) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
- Returns:
Description to print.
- Return type:
str
- str_sym(alat)[source]
Write 4 strings (+1 empty) for the prefix.sym file for sitesym in wannier90: The symmetry operations act on a point r as rR − t.
- Parameters:
alat (float) – Lattice parameter in angstroms.
- Returns:
Description to print. 1 blank line 3 lines: cartesian rotation matrix 1 line : cartesian translation in units of alat
- Return type:
str
- translation_refUC(refUC, shiftUC)[source]
Calculate translation in reference choice of unit cell.
- Parameters:
refUC (array) – 3x3 array describing the transformation of vectors defining the unit cell to the standard setting.
shiftUC (array) – Translation taking the origin of the unit cell used in the DFT calculation to that of the standard setting.
- Returns:
Translation in reference choice of unit cell.
- Return type:
array
- irrep.spacegroup.cart_to_crystal(rot_cart, trans_cart, lattice, alat)[source]
Convert rotation and translation matrices from cartesian to crystal coordinates.
- Parameters:
rot_cart (array, shape=(Nsym, 3, 3)) – Each element is a 3x3 array describing a rotation matrix in cartesian coordinates
trans_cart (array, shape=(Nsym, 3)) – Each element is a 3D vector describing a translation in cartesian coordinates
lattice (array, shape=(3, 3)) – Each row contains cartesian coordinates of a basis vector forming the unit-cell in real space.
alat (float, default=1) – Lattice parameter in angstroms (quantum espresso convention).
- Returns:
array, shape=(Nsym, 3, 3) – Each element is a 3x3 array describing a rotation matrix in crystal coordinates (should be integers)
array, shape=(Nsym, 3) – Each element is a 3D vector describing a translation in crystal coordinates
- irrep.spacegroup.read_sym_file(fname)[source]
Read symmetry operations from a file.
- Parameters:
fname (str) – Name of the file.
- Returns:
np.array(Nsym, 3, 3) – Each element is a 3x3 array describing a rotation matrix in cartesian coordinates
np.array(Nsym, 3) – Each element is a 3D vector describing a translation in units of alat.