Spglib for Python#

Installation#

The source code is downloaded at https://github.com/spglib/spglib/tags. But minor updates are not included in this package. If you want the latest version, you can git-clone the spglib repository

% git clone https://github.com/spglib/spglib.git

It is also possible to install spglib for python via the following package distribution services or from building using setup.py.

The following pip and conda packages are made and maintained by Paweł T. Jochym, which is of great help to keeping spglib handy and useful.

Using pip#

Numpy is required before the python-spglib installation. The command to install spglib is

% pip install --user spglib

If you see the error message like below in the installation process

_spglib.c:35:20: fatal error: Python.h: No such file or directory

development tools for building python module are additionally necessary and are installed using OS’s package management system, e.g.,:

% sudo apt-get install python-dev

If your installation by pip failed, you may need to upgrade setuptools, e.g., by

% pip install --upgrade --user setuptools

Using conda#

Conda is another choice:

% conda install -c conda-forge spglib

Building using setup.py#

To manually install python-spglib using setup.py, python header files (python-dev), C-compiler (e.g., gcc, clang), and numpy are required before the build. The installation steps are shown as follows:

  1. Go to the python directory

  2. Type the command:

    % python setup.py install --user
    

    Document about where spglib is installed is found at the links below:

If your installation by setup.py failed, you may need to upgrade setuptools, e.g., by

% pip install --upgrade --user setuptools

Test#

The test script test_spglib.py is found in python/test directory. Got to this directory and run this script. It will be like below::

% python test_spglib.py
test_find_primitive (__main__.TestSpglib) ... ok
test_get_symmetry (__main__.TestSpglib) ... ok
test_get_symmetry_dataset (__main__.TestSpglib) ... ok
test_refine_cell (__main__.TestSpglib) ... ok

----------------------------------------------------------------------
Ran 4 tests in 13.147s

OK

How to import spglib module#

Change in version 1.9.0!

For versions 1.9.x or later:

import spglib

For versions 1.8.x or before:

from pyspglib import spglib

If the version is not sure:

try:
    import spglib as spg
except ImportError:
    from pyspglib import spglib as spg

Version number#

In version 1.8.3 or later, the version number is obtained by spglib.__version__ or get_version.

Example#

Examples are found in examples directory.

Variables#

Crystal structure (cell)#

A crystal structure is given by a tuple. This tuple format is supported at version 1.9.1 or later.

The tuple format is shown as follows. There are three or four elements in the tuple: cell = (lattice, positions, numbers) or cell = (lattice, positions, numbers, magmoms) where magmoms represents magnetic moments on atoms and is optional.

Lattice parameters lattice are given by a 3x3 matrix with floating point values, where \(\mathbf{a}, \mathbf{b}, \mathbf{c}\) are given as rows, which results in the transpose of the definition for C-API (lattice). Fractional atomic positions positions are given by a Nx3 matrix with floating point values, where N is the number of atoms. Numbers to distinguish atomic species numbers are given by a list of N integers. The magnetic moments magmoms can be specified with get_magnetic_symmetry. magmoms is a list of N floating-point values for collinear cases and a list of Nx3 in cartesian coordinates for non-collinear cases.

lattice = [[a_x, a_y, a_z],
            [b_x, b_y, b_z],
            [c_x, c_y, c_z]]
positions = [[a_1, b_1, c_1],
               [a_2, b_2, c_2],
               [a_3, b_3, c_3],
               ...]
numbers = [n_1, n_2, n_3, ...]
magmoms = [m_1, m_2, m_3, ...]  # Works with get_magnetic_symmetry for a collinear case
# magmoms = [[m_1x, m_1y, m_1z], ...]  # For a non-collinear case

For example, the crystal structure (cell) of L1\(_{2}\)-type AlNi\(_{3}\) is:

lattice = [[1.0, 0.0, 0.0],
            [0.0, 1.0, 0.0],
            [0.0, 0.0, 1.0]]
positions = [[0.0, 0.0, 0.0], # Al
            [0.5, 0.5, 0.0], # Ni
            [0.0, 0.5, 0.5], # Ni
            [0.5, 0.0, 0.5]] # Ni
numbers = [1, 2, 2, 2]        # Al, Ni, Ni, Ni

Version 1.9.5 or later: The methods that use the crystal structure will return None when a crystal structure is not properly given.

ASE Atoms-like input is deprecated.#

In the previous versions, ASE Atoms-like input was supported, but it is deprecated. It is recommended to use the above tuple-style input for the future support. DeprecationWarning is issued at version 1.10.0 or later.

The reason to make this feature deprecated is that ASE Atoms class is too huge and users may expect spglib can understand its full feature. However spglib actually collects only the following values from the ASE Atoms-class instance::

   lattice = cell.get_cell()
   positions = cell.get_scaled_positions()
   numbers = cell.get_atomic_numbers()
   magmoms = None

for which the corresponding code is written out of API and it is found at here. Nevertheless the ASE Atoms-like input will be accepted for a while. An alternative Atoms class (atoms.py) that contains minimum set of methods is prepared in the examples directory. get_symmetry with collinear polarizations is not supported when ASE Atoms-class instance.

Symmetry tolerance (symprec, angle_tolerance, mag_symprec)#

Distance tolerance in Cartesian coordinates to find crystal symmetry. For more details, see symprec, angle_tolerance, and mag_symprec.

Methods (version)#

get_version#

New in version 1.8.3

version = get_version()

This returns version number of spglib by tuple with three numbers.

Methods (error)#

get_error_message#

New in version 1.9.5

Be careful. This method is not thread safe, i.e., only safely usable when calling one spglib method per process.

This method is used to see roughly why spglib failed.

error_message = get_error_message()

Method (standardization and finding primitive)#

standardize_cell#

New in version 1.8.x

lattice, scaled_positions, numbers = standardize_cell(bulk, to_primitive=False, no_idealize=False, symprec=1e-5, angle_tolerance=-1.0)

to_primitive=True is used to create the standardized primitive cell, and no_idealize=True disables to idealize lengths and angles of basis vectors and positions of atoms according to crystal symmetry. Now refine_cell and find_primitive are shorthands of this method with combinations of these options. When the search failed, None is returned. About the default choice of the setting, see the documentation of hall_number argument of get_symmetry_dataset. More detailed explanation is shown in the spglib (C-API) document.

find_primitive#

Behaviour changed in version 1.8.x

lattice, scaled_positions, numbers = find_primitive(cell, symprec=1e-5, angle_tolerance=-1.0)

When a primitive cell is found, lattice parameters (a 3x3 numpy array), scaled positions (a numpy array of [number_of_atoms,3]), and atomic numbers (a 1D numpy array) is returned. When the search failed, None is returned.

The detailed control of standardization of unit cell can be done using standardize_cell.

refine_cell#

Behaviour changed in version 1.8.x

lattice, scaled_positions, numbers = refine_cell(cell, symprec=1e-5, angle_tolerance=-1.0)

Standardized crystal structure is obtained as a tuple of lattice (a 3x3 numpy array), atomic scaled positions (a numpy array of [number_of_atoms,3]), and atomic numbers (a 1D numpy array) that are symmetrized following space group type. When the search failed, None is returned. About the default choice of the setting, see the documentation of hall_number argument of get_symmetry_dataset.

The detailed control of standardization of unit cell is achieved using standardize_cell.

Method (space-group dataset access)#

get_symmetry_from_database#

symmetry = get_symmetry_from_database(hall_number)

A set of crystallographic symmetry operations corresponding to hall_number is returned by a dictionary where rotation parts and translation parts are accessed by the keys rotations and translations, respectively. The definition of hall_number is found at Space group type.

When something wrong happened, None is returned.

get_spacegroup_type#

  • New at version 1.9.4

  • hall_number member is added at version 2.0.

spacegroup_type = get_spacegroup_type(hall_number)

This function allows to directly access to the space-group-type database in spglib (spg_database.c). A dictionary is returned. 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 keys of the returned dictionary is as follows:

number
international_short
international_full
international
schoenflies
hall_symbol
hall_number
choice
pointgroup_schoenflies
pointgroup_international
arithmetic_crystal_class_number
arithmetic_crystal_class_symbol

Here spacegroup_type['international_short'] is equivalent to dataset['international'] of get_symmetry_dataset, spacegroup_type['hall_symbol'] is equivalent to dataset['hall'] of get_symmetry_dataset, and spacegroup_type['pointgroup_international'] is equivalent to dataset['pointgroup_symbol'] of get_symmetry_dataset.

When something wrong happened, None is returned.

get_spacegroup_type_from_symmetry#

New at version 2.0

spacegroup_type = get_spacegroup_type_from_symmetry(rotations, translations, lattice=None, symprec=1e-5)

Return space-group type information (see get_spacegroup_type) from symmetry operations. This is the replacement of get_hall_number_from_symmetry (deprecated at v2.0).

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 parameter lattice is used as the distance measure for symprec. If this is not given, the cubic basis vector whose lengths are one is used.

Methods (magnetic symmetry)#

get_magnetic_symmetry#

Experimental: new at version 2.0

Find magnetic symmetry operations from a crystal structure and site tensors

symmetry = get_magnetic_symmetry(
    cell,
    symprec=1e-5,
    angle_tolerance=-1.0,
    mag_symprec=-1.0,
    is_axial=None,
    with_time_reversal=True,
)

Arguments

  • cell: See Crystal structure (cell)

  • symprec, angle_tolerance: See get_symmetry

  • mag_symprec: Tolerance for magnetic symmetry search in the unit of magmoms. If not specified, use the same value as symprec.

  • is_axial Set is_axial=True if magmoms does not change their sign by improper rotations. If not specified, set is_axial=False when magmoms.shape==(num_atoms, ), and set is_axial=True when magmoms.shape==(num_atoms, 3). These default settings correspond to collinear and non-collinear spins.

  • with_time_reversal Set with_time_reversal=True if magmoms change their sign by time-reversal operations.

Returned symmetry is a dictionary with the following keys:

  • rotations : Rotation (matrix) parts of symmetry operations

  • translations: Translation (vector) parts of symmetry operations

  • time_reversals: Time reversal part of magnetic symmetry operations. True indicates time reversal operation, and False indicates an ordinary operation.

  • equivalent_atoms

get_magnetic_symmetry_dataset#

Experimental: new at version 2.0

Search magnetic symmetry dataset from an input cell.

dataset = get_magnetic_symmetry_dataset(
    cell,
    is_axial=None,
    symprec=1e-5,
    angle_tolerance=-1.0,
    mag_symprec=-1.0,
)

Arguments

Returned dataset is a dictionary. The description of its keys is given at Magnetic Spglib dataset (Experimental).

get_magnetic_spacegroup_type#

Experimental: new at version 2.0

Translate UNI number to magnetic space group type information.

magnetic_spacegroup_type = get_magnetic_spacegroup_type(uni_number):

The returned magnetic_spacegroup_type is a dictionary with

uni_number
litvin_number
bns_number
og_number
number
type

See spg_get_magnetic_spacegroup_type for these descriptions.

get_magnetic_symmetry_from_database#

Experimental: new at version 2.0

symmetry = get_magnetic_symmetry_from_database(uni_number, hall_number=0)

Return magnetic symmetry operations corresponding to uni_number between 1 and 1651. Optionally alternative settings can be specified with hall_number.

Methods (lattice reduction)#

niggli_reduce#

New at version 1.9.4

niggli_lattice = niggli_reduce(lattice, eps=1e-5)

Niggli reduction is achieved using this method. The algorithm detail is found at https://atztogo.github.io/niggli/ and the references are there in. Original basis vectors are stored in lattice and the Niggli reduced basis vectors are given in niggli_lattice. The format of basis vectors are found at Crystal structure (cell). esp is the tolerance parameter, but unlike symprec the unit is not a length. This is used to check if difference of norms of two basis vectors is close to zero or not and if two basis vectors are orthogonal by the value of dot product being close to zero or not. The detail is shown at https://atztogo.github.io/niggli/.

When the search failed, None is returned.

The transformation from original basis vectors \(( \mathbf{a} \; \mathbf{b} \; \mathbf{c} )\) to final basis vectors \(( \mathbf{a}' \; \mathbf{b}' \; \mathbf{c}' )\) is achieved by linear combination of basis vectors with integer coefficients without rotating coordinates. Therefore the transformation matrix is obtained by \(\boldsymbol{P} = ( \mathbf{a} \; \mathbf{b} \; \mathbf{c} ) ( \mathbf{a}' \; \mathbf{b}' \; \mathbf{c}' )^{-1}\) and the matrix elements have to be almost integers.

delaunay_reduce#

New at version 1.9.4

delaunay_lattice = delaunay_reduce(lattice, eps=1e-5)

Delaunay reduction is achieved using this method. The algorithm is found in the international tables for crystallography volume A. Original basis vectors are stored in lattice and the Delaunay reduced basis vectors are given in delaunay_lattice, where the format of basis vectors are shown in Crystal structure (cell). esp is the tolerance parameter, but unlike symprec the unit is not a length. This is used as the criterion if volume is close to zero or not and if two basis vectors are orthogonal by the value of dot product being close to zero or not.

When the search failed, None is returned.

The transformation from original basis vectors \(( \mathbf{a} \; \mathbf{b} \; \mathbf{c} )\) to final basis vectors \(( \mathbf{a}' \; \mathbf{b}' \; \mathbf{c}' )\) is achieved by linear combination of basis vectors with integer coefficients without rotating coordinates. Therefore the transformation matrix is obtained by \(\boldsymbol{P} = ( \mathbf{a} \; \mathbf{b} \; \mathbf{c} ) ( \mathbf{a}' \; \mathbf{b}' \; \mathbf{c}' )^{-1}\) and the matrix elements have to be almost integers.

Methods (kpoints)#

get_ir_reciprocal_mesh#

mapping, grid = get_ir_reciprocal_mesh(mesh, cell, is_shift=[0, 0, 0])

Irreducible k-points are obtained from a sampling mesh of k-points. mesh is given by three integers by array and specifies mesh numbers along reciprocal primitive axis. is_shift is given by the three integers by array. 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. When the value is not 0, is_shift is set.

mapping and grid are returned. grid gives the mesh points in fractional coordinates in reciprocal space. mapping gives mapping to the irreducible k-point indices that are obtained by

np.unique(mapping)

Here np means the numpy module. The grid point is accessed by grid[index].

When the search failed, None is returned.

An example is shown below:

import numpy as np
import spglib

lattice = np.array([[0.0, 0.5, 0.5],
                    [0.5, 0.0, 0.5],
                    [0.5, 0.5, 0.0]]) * 5.4
positions = [[0.875, 0.875, 0.875],
            [0.125, 0.125, 0.125]]
numbers= [1,] * 2
cell = (lattice, positions, numbers)
print(spglib.get_spacegroup(cell, symprec=1e-5))
mesh = [8, 8, 8]

#
# Gamma centre mesh
#
mapping, grid = spglib.get_ir_reciprocal_mesh(mesh, cell, is_shift=[0, 0, 0])

# All k-points and mapping to ir-grid points
for i, (ir_gp_id, gp) in enumerate(zip(mapping, grid)):
    print("%3d ->%3d %s" % (i, ir_gp_id, gp.astype(float) / mesh))

# Irreducible k-points
print("Number of ir-kpoints: %d" % len(np.unique(mapping)))
print(grid[np.unique(mapping)] / np.array(mesh, dtype=float))

#
# With shift
#
mapping, grid = spglib.get_ir_reciprocal_mesh(mesh, cell, is_shift=[1, 1, 1])

# All k-points and mapping to ir-grid points
for i, (ir_gp_id, gp) in enumerate(zip(mapping, grid)):
    print("%3d ->%3d %s" % (i, ir_gp_id, (gp + [0.5, 0.5, 0.5]) / mesh))

# Irreducible k-points
print("Number of ir-kpoints: %d" % len(np.unique(mapping)))
print((grid[np.unique(mapping)] + [0.5, 0.5, 0.5]) / mesh)

Method (deprecated)#

get_hall_number_from_symmetry#

Deprecated at version 2.0. This function is replaced by get_spacegroup_type_from_symmetry.

Return one of hall_number corresponding to a space-group type of the given set of symmetry operations. When multiple hall_number exist for the space-group type, the smallest one (the first description of the space-group type in International Tables for Crystallography) is chosen. The definition of hall_number is found at Space group type and the corresponding space-group-type information is obtained through get_spacegroup_type.

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.

get_hall_number_from_symmetry(rotations, translations, symprec=1e-5)