Public API

Contents

Index

Public interface

To utilize the public interface, first import the required packages:

julia> using CrystallographyCore, Spglib

For documentation on CrystallographyCore.jl, refer to this link.

For extended functionalities, consider exploring my other packages: CrystallographyBase.jl, Crystallography.jl, and MillerIndices.jl.

Types

Spglib.DatasetType
Dataset(spacegroup_number, hall_number, international_symbol, hall_symbol, choice, transformation_matrix, origin_shift, n_operations, rotations, translations, n_atoms, wyckoffs, site_symmetry_symbols, equivalent_atoms, crystallographic_orbits, primitive_lattice, mapping_to_primitive, n_std_atoms, std_lattice, std_types, std_positions, std_rotation_matrix, std_mapping_to_primitive, pointgroup_symbol)

Represent SpglibDataset, see its official documentation.

Arguments

  • spacegroup_number: international space group number.
  • hall_number: Hall number. This number is used in get_symmetry_from_database and get_spacegroup_type.
  • international_symbol: international short symbol.
  • hall_symbol: Hall symbol.
  • choice: centring, origin, basis vector setting.
  • transformation_matrix: See the detail at Transformation matrix and origin shift.
  • origin shift: See the detail at Transformation matrix and origin shift.
  • n_operations: number of symmetry operations.
  • rotations and translations: rotation matrices and translation vectors. See get_symmetry for more details.
  • n_atoms: number of atoms in the input unit cell.
  • wyckoffs: Wyckoff letters.
  • site_symmetry_symbols: site-symmetry symbols (experimental).
  • equivalent_atoms: mapping table to equivalent atoms.
  • crystallographic_orbits : mapping table to equivalent atoms (see Wyckoff positions and symmetrically equivalent atoms for the difference between equivalent_atoms and crystallographic_orbits).
  • primitive_lattice : basis vectors of a primitive cell.
  • mapping_to_primitive: mapping table to atoms in the primitive cell.
  • n_std_atoms: number of atoms in the standardized unit cell.
  • std_lattice, std_positions, std_types: standardized crystal structure corresponding to the Hall symbol found. These are equivalently given in the array formats of lattice, positions, and atoms presented at SpglibCell, respectively.
  • std_rotation_matrix: see the detail at Standardized crystal structure after idealization.
  • std_mapping_to_primitive: Mapping table from atoms in the standardized crystal structure to the atoms in the primitive cell.
  • pointgroup_symbol: symbol of the crystallographic point group in the Hermann–Mauguin notation.

See also get_dataset, get_dataset_with_hall_number.

source
Spglib.MagneticDatasetType
MagneticDataset(uni_number, msg_type, hall_number, tensor_rank, n_operations, rotations, translations, time_reversals, n_atoms, equivalent_atoms, transformation_matrix, origin_shift, n_std_atoms, std_lattice, std_types, std_positions, std_tensors, std_rotation_matrix, primitive_lattice)

Represent MagneticDataset, see its official documentation.

See also get_magnetic_dataset.

source
Spglib.MagneticSpacegroupTypeType

MagneticSpacegroupType(uninumber, litvinnumber, bnsnumber, ognumber, number, type)

Represent SpglibMagneticSpacegroupType, see its official documentation.

Arguments

  • uni_number::Int32: Serial number of UNI (or BNS) symbols.
  • litvin_number::Int32: Serial number in Litvin's Magnetic Group Tables.
  • bns_number::String: BNS number, e.g. "151.32".
  • og_number::String: OG number, e.g. "153.4.1270".
  • number::Int32: ITA's serial number of space group for reference setting.
  • type::Int32: Type of MSG, from $1$ to $4$.

See also get_magnetic_spacegroup_type, get_magnetic_spacegroup_type_from_symmetry.

source
Spglib.SpglibCellType
SpglibCell(lattice, positions, atoms, magmoms=[])

Represent a unit cell with specified lattice, positions, atoms, and magnetic moments.

Arguments

  • lattice: lattice of the unit cell. Lattice parameters are given by a $3×3$ matrix with floating point values, where $𝐚$, $𝐛$, and $𝐜$ are stored as columns. You could also give a vector of 3-vectors, where each vector is a lattice vector. See Basis vectors for our conventions and Lattice for more examples.
  • positions: positions of the atoms in the unit cell. Fractional atomic positions are given by a vector of $N$ 3-vectors with floating point values, where $N$ is the number of atoms.
  • atoms: $N$ atoms present in the unit cell.
  • magmoms=[]: magnetic moments on atoms in the unit cell (optional). It can be either a vector of $N$ floating point values for collinear cases or a vector of 3-vectors in cartesian coordinates for non-collinear cases.

See also Lattice.

Examples

julia> lattice = Lattice([
           [5.0759761474456697, 5.0759761474456697, 0],
           [-2.8280307701821314, 2.8280307701821314, 0],
           [0, 0, 8.57154746],
       ]);

julia> positions = [
           [0.0, 0.84688439, 0.1203133],
           [0.0, 0.65311561, 0.6203133],
           [0.0, 0.34688439, 0.3796867],
           [0.0, 0.15311561, 0.8796867],
           [0.5, 0.34688439, 0.1203133],
           [0.5, 0.15311561, 0.6203133],
           [0.5, 0.84688439, 0.3796867],
           [0.5, 0.65311561, 0.8796867],
       ];

julia> atoms = fill(35, length(positions));

julia> cell = SpglibCell(lattice, positions, atoms);

julia> lattice = Lattice([
           4 0 0
           0 4 0
           0 0 3
       ]);

julia> positions = [
           [0.0, 0.0, 0.0],
           [0.5, 0.5, 0.5],
           [0.3, 0.3, 0.0],
           [0.7, 0.7, 0.0],
           [0.2, 0.8, 0.5],
           [0.8, 0.2, 0.5],
       ];

julia> atoms = [14, 14, 8, 8, 8, 8];

julia> cell = SpglibCell(lattice, positions, atoms);

julia> lattice = [
           4.0 0.0 0.0
           0.0 4.0 0.0
           0.0 0.0 4.0
       ];

julia> positions = [[0.0, 0.0, 0.0], [0.5, 0.5, 0.5]];

julia> atoms = [1, 1];

julia> magmoms = [1.0, 1.0];

julia> cell = SpglibCell(lattice, positions, atoms, magmoms);
source

Functions

You can find their official documentation on this page.

Spglib.delaunay_reduceFunction
delaunay_reduce(lattice::AbstractMatrix, symprec=1e-5)
delaunay_reduce(cell::Cell, symprec=1e-5)

Apply Delaunay reduction to input basis vectors lattice.

The transformation from original basis vectors $\begin{bmatrix} \mathbf{a} & \mathbf{b} & \mathbf{c} \end{bmatrix}$ to final basis vectors $\begin{bmatrix} \mathbf{a}' & \mathbf{b}' & \mathbf{c}' \end{bmatrix}$ is achieved by linear combination of basis vectors with integer coefficients without rotating coordinates. Therefore the transformation matrix is obtained by

\[\mathbf{P} = \begin{bmatrix} \mathbf{a} & \mathbf{b} & \mathbf{c} \end{bmatrix} \begin{bmatrix} \mathbf{a}' & \mathbf{b}' & \mathbf{c}' \end{bmatrix}^{-1}\]

and the matrix elements have to be almost integers.

See also Computing rigid rotation introduced by idealization.

source
Spglib.eachpointFunction
eachpoint(result::BrillouinZoneMesh, ir_only=true)

Iterate over the points in the Brillouin zone mesh, with the option to include only irreducible k-points or all k-points.

See also BrillouinZoneMesh.

Note

This function only returns an iterator, not a vector of points. To get a vector, use collect(eachpoint(result, ir_only)).

source
Spglib.find_primitiveFunction
find_primitive(cell::AbstractCell, symprec=1e-5)

Find the primitive cell of an input unit cell.

This function is now a shortcut of standardize_cell with to_primitive=true and no_idealize=false.

source
Spglib.get_datasetFunction
get_dataset(cell::AbstractCell, symprec=1e-5)

Search symmetry operations of an input unit cell structure.

For an input unit cell structure, symmetry operations of the crystal are searched. Then they are compared with the crystallographic database and the space group type is determined. The result is returned as the Dataset structure as a dataset.

The detail of the dataset is given at Dataset.

Dataset corresponding to the space group type in the standard setting is obtained by get_dataset. Here the standard setting means the first top one among the Hall symbols listed for each space group type. For example, H setting (hexagonal lattice) is chosen for rhombohedral crystals. get_dataset_with_hall_number explained below is used to choose different settings such as R setting of rhombohedral crystals. In this function, the other crystallographic setting is not obtained.

source
Spglib.get_dataset_with_hall_numberFunction
get_dataset_with_hall_number(cell::AbstractCell, hall_number::Integer, symprec=1e-5)

Search symmetry operations of an input unit cell structure, using a given Hall number.

For an input unit cell structure, symmetry operations of the crystal are searched. Then they are compared with the crystallographic database and the space group type is determined. The result is returned as the Dataset structure as a dataset.

The detail of the dataset is given at Dataset.

Dataset corresponding to the space group type in the standard setting is obtained by get_dataset. Here the standard setting means the first top one among the Hall symbols listed for each space group type. For example, H setting (hexagonal lattice) is chosen for rhombohedral crystals. get_dataset_with_hall_number explained below is used to choose different settings such as R setting of rhombohedral crystals. In this function, the other crystallographic setting is not obtained.

source
Spglib.get_hall_number_from_symmetryFunction
get_hall_number_from_symmetry(rotations, translations, symprec=1e-5)

Return one Hall number corresponding to a space group of the given set of symmetry operations.

When multiple Hall numbers exist for the space group, the smallest one (the first description of the space-group type in International Tables for Crystallography) is chosen.

This is expected to work well for the set of symmetry operations whose distortion is small. The aim of making this feature is to find space-group-type for the set of symmetry operations given by the other source than Spglib.

Note that the definition of symprec is different from usual one, but is given in the fractional coordinates and so it should be small like 1e-5.

Warning

This function will be replaced by get_spacegroup_type_from_symmetry.

source
Spglib.get_internationalFunction
get_international(cell::AbstractCell, symprec=1e-5)

Return the space group type in Hermann–Mauguin (international) notation.

source
Spglib.get_ir_reciprocal_meshFunction
get_ir_reciprocal_mesh(cell::AbstractCell, mesh, symprec=1e-5; is_shift=falses(3), is_time_reversal=true)

Search irreducible reciprocal grid points from uniform mesh grid points specified by mesh and is_shift.

Reciprocal primitive vectors are divided by the number stored in mesh with $(0, 0, 0)$-centering. The center of grid mesh is shifted half of a grid spacing along corresponding reciprocal axis by setting 1 or true to each is_shift element. If 0 or false is set to each is_shift element, there is no shift. This limitation of shifting enables the irreducible k-point search significantly faster when the mesh is very dense.

Arguments

  • cell: the input cell.
  • mesh: the mesh numbers along each reciprocal axis. It is given by three integers.
  • symprec: the tolerance for symmetry search.
  • is_shift: a 3-Boolean vector. When is_shift is set for each reciprocal primitive axis, the mesh is shifted along the axis in half of adjacent mesh points irrespective of the mesh numbers.
  • is_time_reversal: whether to impose the time reversal symmetry.

See also BrillouinZoneMesh.

source
Spglib.get_magnetic_spacegroup_type_from_symmetryFunction
get_magnetic_spacegroup_type_from_symmetry(rotations, translations, time_reversals, lattice::Lattice, symprec=1e-5)

Determine magnetic space-group type from magnetic symmetry operations.

time_reversals takes false for ordinary operations and true for time-reversal operations.

source
Spglib.get_magnetic_symmetry_from_databaseFunction
get_magnetic_symmetry_from_database(uni_number, hall_number=0)

Accesses magnetic space-group operations in the built-in database using UNI number (from $1$ to $1651$).

Optionally alternative settings can be specified with hallnumber. For type-I, type-II, and type-III magnetic space groups, `hallnumberchanges settings in family space group. For type-IV,hallnumberchanges settings in maximal space group. Whenhallnumber=0, the smallest hall number corresponding touni_number` is used.

source
Spglib.get_multiplicityFunction
get_multiplicity(cell::AbstractCell, symprec=1e-5)

Return the exact number of symmetry operations.

An error is thrown when it fails.

source
Spglib.get_schoenfliesFunction
get_schoenflies(cell::AbstractCell, symprec=1e-5)

Return the space group type in Schoenflies notation.

source
Spglib.get_spacegroup_type_from_symmetryFunction
get_spacegroup_type_from_symmetry(rotations, translations, lattice::Lattice, symprec=1e-5)

Return space-group type information from symmetry operations.

This is the replacement of get_hall_number_from_symmetry.

This is expected to work well for the set of symmetry operations whose distortion is small. The aim of making this feature is to find space-group-type for the set of symmetry operations given by the other source than Spglib.

The SpacegroupType structure is explained at SpacegroupType. The parameter lattice is used as the distance measure for symprec. If it is unknown, the following may be a reasonable choice:

julia> lattice = Lattice([
           1 0 0
           0 1 0
           0 0 1
       ]);
source
Spglib.get_stabilized_reciprocal_meshFunction
get_stabilized_reciprocal_mesh(rotations, mesh, qpoints=[[0, 0, 0]]; is_shift=falses(3), is_time_reversal=true)

Search irreducible k-points from unique k-point mesh grids from direct (real space) basis vectors and a set of rotation parts of symmetry operations in direct space with one or multiple stabilizers.

The stabilizers are written in fractional coordinates. This function can be used to obtain all mesh grid points by setting rotations = [[1 0 0; 0 1 0; 0 0 1]], and qpoints = [[0, 0, 0]].

See also BrillouinZoneMesh.

source
Spglib.get_symmetryFunction
get_symmetry(cell::AbstractCell, symprec=1e-5)

Return the symmetry operations (rotations, translations) of a cell.

Returned value rotations is a Vector of matrices. It has the length of "number of symmetry operations". Each matrix is a $3 \times 3$ integer matrix. Returned value translations is a Vector of vectors. It has the length of "number of symmetry operations". Each vector is a length-$3$ vector of floating point numbers.

The orders of the rotation matrices and the translation vectors correspond with each other, e.g., the second symmetry operation is organized by the set of the second rotation matrix and second translation vector in the respective arrays. Therefore a set of symmetry operations may obtained by [(r, t) for r, t in zip(rotations, translations)].

The operations are given with respect to the fractional coordinates (not for Cartesian coordinates). The rotation matrix $\mathbf{W}$ and translation vector $\text{w}$ are used as follows:

\[\tilde{\mathbf{x}}_{3\times 1} = \mathbf{W}_{3\times 3} \mathbf{x}_{3\times 1} + \text{w}_{3\times 1}.\]

The three values in the vector are given for the $a$, $b$, and $c$ axes, respectively.

As an exceptional case, if a supercell (or non-primitive cell) has the basis vectors whose lattice breaks crystallographic point group, the crystallographic symmetry operations are searched within this broken symmetry, i.e., at most the crystallographic point group found in this case is the point group of the lattice. For example, this happens for the $2\times 1\times 1$ supercell of a conventional cubic unit cell. This may not be understandable in crystallographic sense, but is practically useful treatment for research in computational materials science.

Examples

julia> lattice = Lattice([
           [5.0759761474456697, 5.0759761474456697, 0],
           [-2.8280307701821314, 2.8280307701821314, 0],
           [0, 0, 8.57154746],
       ]);

julia> positions = [
           [0.0, 0.84688439, 0.1203133],
           [0.0, 0.65311561, 0.6203133],
           [0.0, 0.34688439, 0.3796867],
           [0.0, 0.15311561, 0.8796867],
           [0.5, 0.34688439, 0.1203133],
           [0.5, 0.15311561, 0.6203133],
           [0.5, 0.84688439, 0.3796867],
           [0.5, 0.65311561, 0.8796867],
       ];

julia> atoms = fill(35, length(positions));

julia> cell = SpglibCell(lattice, positions, atoms);

julia> rotations, translations = get_symmetry(cell);

julia> size(rotations) == size(translations) == (16,)
true
source
Spglib.get_symmetry_from_databaseMethod
get_symmetry_from_database(hall_number)

Return the symmetry operations given a hall_number.

This function allows to directly access to the space group operations in the spglib database. To specify the space group type with a specific choice, hall_number is used.

The definition of hall_number is found at Space group type.

The serial number from $1$ to $530$ which are found at list of space groups (Seto's web site). Be sure that this is not a standard crystallographic definition.

Examples

julia> get_symmetry_from_database(304);
source
Spglib.get_symmetry_with_collinear_spinFunction
get_symmetry_with_collinear_spin(cell::SpglibCell, symprec=1e-5)

Find symmetry operations with collinear polarizations (spins) on atoms.

Except for the magmoms in the cell, the usage is basically the same as get_symmetry. But as an output, equivalent_atoms are obtained as the last returned value. The size of this array is the same of number of atoms in the cell.

source
Spglib.get_symmetry_with_site_tensorsFunction
get_symmetry_with_site_tensors(cell::SpglibCell, symprec=1e-5; with_time_reversal=true)

Return magnetic symmetry operations represented by rotation, translation, and spin_flips.

Returned spin_flips represents sign of site tensors after applying time-reversal operations: $1$ for non time reversal, and $-1$ for time reversal.

source
Spglib.get_versionMethod
get_version()

Obtain the version number of spglib.

This is the mergence of spg_get_major_version, spg_get_minor_version, and spg_get_micro_version in its C-API.

source
Spglib.niggli_reduceFunction
niggli_reduce(lattice::AbstractMatrix, symprec=1e-5)
niggli_reduce(cell::Cell, symprec=1e-5)

Apply Niggli reduction to input basis vectors lattice.

The transformation from original basis vectors $\begin{bmatrix} \mathbf{a} & \mathbf{b} & \mathbf{c} \end{bmatrix}$ to final basis vectors $\begin{bmatrix} \mathbf{a}' & \mathbf{b}' & \mathbf{c}' \end{bmatrix}$ is achieved by linear combination of basis vectors with integer coefficients without rotating coordinates. Therefore the transformation matrix is obtained by

\[\mathbf{P} = \begin{bmatrix} \mathbf{a} & \mathbf{b} & \mathbf{c} \end{bmatrix} \begin{bmatrix} \mathbf{a}' & \mathbf{b}' & \mathbf{c}' \end{bmatrix}^{-1}\]

and the matrix elements have to be almost integers.

See also Computing rigid rotation introduced by idealization.

source
Spglib.refine_cellFunction
refine_cell(cell::AbstractCell, symprec=1e-5)

Return the refined cell.

This function is now a shortcut of standardize_cell with to_primitive=false and no_idealize=false.

The standardized crystal structure is obtained from a non-standard crystal structure which may be slightly distorted within a symmetry recognition tolerance, or whose primitive vectors are differently chosen, etc.

source
Spglib.standardize_cellFunction
standardize_cell(cell::AbstractCell, symprec=1e-5; to_primitive=false, no_idealize=false)

Return the standardized cell.

The standardized unit cell (see Spglib conventions of standardized unit cell) is generated from an input unit cell structure and its symmetry found by the symmetry search. The choice of the setting for each space group type is as explained for get_dataset. Usually to_primitive=false and no_idealize=false are recommended to set and this setting results in the same behavior as spg_refine_cell.

The standardized unit cell is generated from an input unit cell structure and its symmetry found by the symmetry search. The choice of the setting for each space group type is as explained for get_dataset.

Arguments

  • cell: the input cell to standardize.
  • symprec: the tolerance for symmetry search.
  • to_primitive=true is used to create the standardized primitive cell with the transformation matrices shown at Transformation to the primitive cell, otherwise to_primitive=false must be specified.
  • no_idealize=false is used to idealize the lengths and angles of basis vectors with adjusting the positions of atoms to nearest exact positions according to crystal symmetry. However the crystal can be rotated in Cartesian coordinates by the idealization of the basis vectors. no_idealize=true disables this. The detail of the idealization (no_idealize=false) is written at Idealization of unit cell structure. no_idealize=true may be useful when we want to leave basis vectors and atomic positions in Cartesian coordinates fixed.
source