Troubleshooting

This page collects some possible errors you may encounter along with tips on how to fix them. If you have some questions about how to use this code, you are welcome to discuss with us.

If you have additional tips, please either report an issue or submit a pull request with suggestions.

Cannot find the Julia executable

Make sure you have Julia installed in your environment. Please download the latest stable version for your platform. If you are using a *nix system, the recommended way is to use Juliaup. If you do not want to install Juliaup or you are using other platforms that Julia supports, download the corresponding binaries. Then, create a symbolic link to the Julia executable. If the path is not in your $PATH environment variable, export it to your $PATH.

Some clusters, like Comet, or Expanse, already have Julia installed as a module, you may just module load julia to use it. If not, either install by yourself or contact your administrator.

See Installation Guide for more information.

Julia starts slow

First, we recommend you download the latest version of Julia. Usually, the newest version has the best performance.

If you need to use Julia for a simple, one-time task, you can start the Julia REPL with

julia --compile=min

to minimize compilation or

julia --optimize=0

to minimize optimizations, or just use both. Or you could make a system image and run with

julia --sysimage custom-image.so

See Fredrik Ekre's talk for details.

Unable to install the package

First, verify if your platform is supported by referring to the spglib_jll.jl platforms. If it's not supported, kindly report the issue to the developers here.

If you can install spglib_jll.jl but encounter issues with Spglib.jl, please report to the developer here.

Returned cell symmetry is wrong

Check whether you set the lattice correctly. This is the part where errors can easily occur, as we adopt a different convention from the Python and C versions.

For example, the example shown in Spglib's official documentation is written as follows:

#include <assert.h>
#include "spglib.h"

int main(void) {
    SpglibDataset *dataset;
    // Wurtzite structure (P6_3mc)
    double lattice[3][3] = {
        {3.111, -1.5555, 0}, {0, 2.6942050311733885, 0}, {0, 0, 4.988}};
    double position[4][3] = {
        {1.0 / 3, 2.0 / 3, 0.0},
        {2.0 / 3, 1.0 / 3, 0.5},
        {1.0 / 3, 2.0 / 3, 0.6181},
        {2.0 / 3, 1.0 / 3, 0.1181},
    };
    int types[4] = {1, 1, 2, 2};
    int num_atom = 4;
    double symprec = 1e-5;
    dataset = spg_get_dataset(lattice, position, types, num_atom, symprec);
    assert(dataset->spacegroup_number == 186);
    spg_free_dataset(dataset);
}

Thus, the Python correspondence of the code should be:

import spglib

lattice = [[3.111, 0, 0], [-1.5555, 2.6942050311733885, 0], [0, 0, 4.988]]
position = [
    [1.0 / 3, 2.0 / 3, 0.0],
    [2.0 / 3, 1.0 / 3, 0.5],
    [1.0 / 3, 2.0 / 3, 0.6181],
    [2.0 / 3, 1.0 / 3, 0.1181],
]
types = [1, 1, 2, 2]
symprec = 1e-5
cell = (lattice, positions, types)
dataset = spglib.get_symmetry_dataset(cell, symprec=symprec)
assert dataset['number'] == 186

Note that in Python, the lattice is transposed, as explained in its official documentation.

However, the corresponding code in Julia should be written as follows:

julia> using Spglib
julia> lattice = [[3.111, 0, 0], [-1.5555, 2.6942050311733885, 0], [0, 0, 4.988]];
julia> positions = [ [1.0 / 3, 2.0 / 3, 0.0], [2.0 / 3, 1.0 / 3, 0.5], [1.0 / 3, 2.0 / 3, 0.6181], [2.0 / 3, 1.0 / 3, 0.1181], ];
julia> atoms = [1, 1, 2, 2];
julia> cell = Cell(lattice, positions, atoms)SpglibCell{Float64, Float64, Int64, Any} lattice: 3.111 -1.5555 0.0 0.0 2.6942050311733885 0.0 0.0 0.0 4.988 4 atomic positions: 0.3333333333333333 0.6666666666666666 0.0 0.6666666666666666 0.3333333333333333 0.5 0.3333333333333333 0.6666666666666666 0.6181 0.6666666666666666 0.3333333333333333 0.1181 4 atoms: 1 1 2 2
julia> dataset = get_dataset(cell, 1e-5)Dataset spacegroup_number: 186 hall_number: 480 international_symbol: P6_3mc hall_symbol: P 6c -2c choice: transformation_matrix: [1.0 -5.551115123125783e-17 0.0; 0.0 0.9999999999999999 0.0; 0.0 0.0 0.9999999999999999] origin_shift: [-5.551115123125783e-17, 0.0, 0.0] n_operations: 12 rotations: StaticArraysCore.SMatrix{3, 3, Int32, 9}[[1 0 0; 0 1 0; 0 0 1], [1 -1 0; 1 0 0; 0 0 1], [0 -1 0; 1 -1 0; 0 0 1], [-1 0 0; 0 -1 0; 0 0 1], [-1 1 0; -1 0 0; 0 0 1], [0 1 0; -1 1 0; 0 0 1], [0 1 0; 1 0 0; 0 0 1], [1 0 0; 1 -1 0; 0 0 1], [1 -1 0; 0 -1 0; 0 0 1], [0 -1 0; -1 0 0; 0 0 1], [-1 0 0; -1 1 0; 0 0 1], [-1 1 0; 0 1 0; 0 0 1]] translations: StaticArraysCore.SVector{3, Float64}[[0.0, 0.0, 0.0], [3.0814879110195774e-33, -5.551115123125782e-17, 0.5], [5.551115123125783e-17, -5.551115123125782e-17, 0.0], [1.1102230246251565e-16, 0.0, 0.5], [1.1102230246251565e-16, 5.551115123125782e-17, 0.0], [5.551115123125783e-17, 5.551115123125782e-17, 0.5], [5.551115123125783e-17, -5.551115123125782e-17, 0.5], [3.0814879110195774e-33, -5.551115123125782e-17, 0.0], [0.0, 0.0, 0.5], [5.551115123125783e-17, 5.551115123125782e-17, 0.0], [1.1102230246251565e-16, 5.551115123125782e-17, 0.5], [1.1102230246251565e-16, 0.0, 0.0]] n_atoms: 4 wyckoffs: ['b', 'b', 'b', 'b'] site_symmetry_symbols: ["3m.", "3m.", "3m.", "3m."] equivalent_atoms: Int32[1, 1, 3, 3] crystallographic_orbits: Int32[1, 1, 3, 3] primitive_lattice: [3.111 -1.5555 0.0; 0.0 2.6942050311733885 0.0; 0.0 0.0 4.988] mapping_to_primitive: Int32[1, 2, 3, 4] n_std_atoms: 4 std_lattice: [3.111 -1.5555 0.0; 0.0 2.6942050311733885 0.0; 0.0 0.0 4.988] std_types: Int32[1, 1, 2, 2] std_positions: StaticArraysCore.SVector{3, Float64}[[0.3333333333333333, 0.6666666666666666, 0.0], [0.6666666666666667, 0.3333333333333333, 0.5], [0.3333333333333333, 0.6666666666666666, 0.6181], [0.6666666666666667, 0.3333333333333333, 0.1181000000000001]] std_rotation_matrix: [1.0 0.0 0.0; 0.0 1.0 0.0; 0.0 0.0 1.0] std_mapping_to_primitive: Int32[1, 2, 3, 4] pointgroup_symbol: 6mm
julia> dataset.spacegroup_number186

Although the Julia definition of our Basis vectors is not transposed (like the C-API), when written one by one, it still resembles a transposed version of the lattice in C. However, they both represent the following matrix:

\[\begin{bmatrix} 3.111 & -1.5555 & 0 \\ 0 & 2.6942050311733885 & 0 \\ 0 & 0 & 4.988 \end{bmatrix}\]

Of course, you can construct the lattice directly using its matrix form:

julia> lattice = [
           3.111  -1.5555  0.0
           0.0  2.6942050311733885  0.0
           0.0  0.0  4.988
       ];

If you are not careful when writing these matrices, you may encounter unexpected results:

import spglib

lattice = [
    [3.111, -1.5555, 0],
    [0, 2.6942050311733885, 0],
    [0, 0, 4.988]
]
positions = [
    [1.0 / 3, 2.0 / 3, 0.0],
    [2.0 / 3, 1.0 / 3, 0.5],
    [1.0 / 3, 2.0 / 3, 0.6181],
    [2.0 / 3, 1.0 / 3, 0.1181]
]
types = [1, 1, 2, 2]
num_atom = 4
symprec = 1e-5
cell = (lattice, positions, types)
dataset = spglib.get_symmetry_dataset(cell, symprec=symprec)

>>> dataset['number']
4
>>> dataset['international']
'P2_1'

Or this:

julia> lattice = [[3.111, -1.5555, 0], [0, 2.6942050311733885, 0], [0, 0, 4.988]];
julia> cell = Cell(lattice, positions, atoms)SpglibCell{Float64, Float64, Int64, Any} lattice: 3.111 0.0 0.0 -1.5555 2.6942050311733885 0.0 0.0 0.0 4.988 4 atomic positions: 0.3333333333333333 0.6666666666666666 0.0 0.6666666666666666 0.3333333333333333 0.5 0.3333333333333333 0.6666666666666666 0.6181 0.6666666666666666 0.3333333333333333 0.1181 4 atoms: 1 1 2 2
julia> dataset = get_dataset(cell, 1e-5)Dataset spacegroup_number: 4 hall_number: 6 international_symbol: P2_1 hall_symbol: P 2yb choice: b transformation_matrix: [-0.9999999999999999 0.9999999999999999 0.0; -0.0 0.0 -0.9999999999999999; -1.0 0.0 0.0] origin_shift: [0.0, 0.0, 0.0] n_operations: 2 rotations: StaticArraysCore.SMatrix{3, 3, Int32, 9}[[1 0 0; 0 1 0; 0 0 1], [-1 0 0; 0 -1 0; 0 0 1]] translations: StaticArraysCore.SVector{3, Float64}[[-0.0, -0.0, -0.0], [-0.0, -0.0, 0.5]] n_atoms: 4 wyckoffs: ['a', 'a', 'a', 'a'] site_symmetry_symbols: ["1", "1", "1", "1"] equivalent_atoms: Int32[1, 1, 3, 3] crystallographic_orbits: Int32[1, 1, 3, 3] primitive_lattice: [-0.0 -3.111 -0.0; -2.6942050311733885 -1.1387050311733884 -0.0; -0.0 -0.0 -4.988] mapping_to_primitive: Int32[1, 2, 3, 4] n_std_atoms: 4 std_lattice: [2.6942050311733885 0.0 -1.1387050311733886; 0.0 4.988 0.0; 0.0 0.0 3.1110000000000007] std_types: Int32[1, 1, 2, 2] std_positions: StaticArraysCore.SVector{3, Float64}[[0.3333333333333336, 0.0, 0.6666666666666667], [0.6666666666666664, 0.5, 0.33333333333333326], [0.3333333333333336, 0.3818999999999999, 0.6666666666666667], [0.6666666666666664, 0.8818999999999999, 0.33333333333333326]] std_rotation_matrix: [0.0 1.0 0.0; 0.0 0.0 -1.0; -1.0 0.0 0.0] std_mapping_to_primitive: Int32[1, 2, 3, 4] pointgroup_symbol: 2
julia> dataset.spacegroup_number4
julia> dataset.international_symbol"P2_1"

Julia's results are different from Python's

For the same reason, the returned results from Python, especially lattices, are transposed versions of those in Julia.

>>> dataset['primitive_lattice']
array([[ 3.111     ,  0.        ,  0.        ],
       [-1.5555    ,  2.69420503,  0.        ],
       [ 0.        ,  0.        ,  4.988     ]])
>>> dataset['std_lattice']
array([[ 3.111     ,  0.        ,  0.        ],
       [-1.5555    ,  2.69420503,  0.        ],
       [ 0.        ,  0.        ,  4.988     ]])
>>> std_lattice_before_idealization = np.dot(
        np.transpose(lattice),
        np.linalg.inv(dataset['transformation_matrix'])
    ).T
array([[ 3.111     ,  0.        ,  0.        ],
       [-1.5555    ,  2.69420503,  0.        ],
       [ 0.        ,  0.        ,  4.988     ]])

These are equivalent to the following Julia matrices:

julia> [
           3.111   0.0      0.0
           -1.5555  2.69421  0.0
           0.0     0.0      4.988
       ];
julia> [ 3.111 0.0 0.0 -1.5555 2.69421 0.0 0.0 0.0 4.988 ];
julia> [ 3.111 0.0 0.0 -1.5555 2.69421 0.0 0.0 0.0 4.988 ];

However, the actual Julia results should be:

julia> using Spglib
julia> lattice = [[3.111, 0, 0], [-1.5555, 2.6942050311733885, 0], [0, 0, 4.988]];
julia> positions = [ [1.0 / 3, 2.0 / 3, 0.0], [2.0 / 3, 1.0 / 3, 0.5], [1.0 / 3, 2.0 / 3, 0.6181], [2.0 / 3, 1.0 / 3, 0.1181], ];
julia> atoms = [1, 1, 2, 2];
julia> cell = Cell(lattice, positions, atoms)SpglibCell{Float64, Float64, Int64, Any} lattice: 3.111 -1.5555 0.0 0.0 2.6942050311733885 0.0 0.0 0.0 4.988 4 atomic positions: 0.3333333333333333 0.6666666666666666 0.0 0.6666666666666666 0.3333333333333333 0.5 0.3333333333333333 0.6666666666666666 0.6181 0.6666666666666666 0.3333333333333333 0.1181 4 atoms: 1 1 2 2
julia> dataset = get_dataset(cell, 1e-5)Dataset spacegroup_number: 186 hall_number: 480 international_symbol: P6_3mc hall_symbol: P 6c -2c choice: transformation_matrix: [1.0 -5.551115123125783e-17 0.0; 0.0 0.9999999999999999 0.0; 0.0 0.0 0.9999999999999999] origin_shift: [-5.551115123125783e-17, 0.0, 0.0] n_operations: 12 rotations: StaticArraysCore.SMatrix{3, 3, Int32, 9}[[1 0 0; 0 1 0; 0 0 1], [1 -1 0; 1 0 0; 0 0 1], [0 -1 0; 1 -1 0; 0 0 1], [-1 0 0; 0 -1 0; 0 0 1], [-1 1 0; -1 0 0; 0 0 1], [0 1 0; -1 1 0; 0 0 1], [0 1 0; 1 0 0; 0 0 1], [1 0 0; 1 -1 0; 0 0 1], [1 -1 0; 0 -1 0; 0 0 1], [0 -1 0; -1 0 0; 0 0 1], [-1 0 0; -1 1 0; 0 0 1], [-1 1 0; 0 1 0; 0 0 1]] translations: StaticArraysCore.SVector{3, Float64}[[0.0, 0.0, 0.0], [3.0814879110195774e-33, -5.551115123125782e-17, 0.5], [5.551115123125783e-17, -5.551115123125782e-17, 0.0], [1.1102230246251565e-16, 0.0, 0.5], [1.1102230246251565e-16, 5.551115123125782e-17, 0.0], [5.551115123125783e-17, 5.551115123125782e-17, 0.5], [5.551115123125783e-17, -5.551115123125782e-17, 0.5], [3.0814879110195774e-33, -5.551115123125782e-17, 0.0], [0.0, 0.0, 0.5], [5.551115123125783e-17, 5.551115123125782e-17, 0.0], [1.1102230246251565e-16, 5.551115123125782e-17, 0.5], [1.1102230246251565e-16, 0.0, 0.0]] n_atoms: 4 wyckoffs: ['b', 'b', 'b', 'b'] site_symmetry_symbols: ["3m.", "3m.", "3m.", "3m."] equivalent_atoms: Int32[1, 1, 3, 3] crystallographic_orbits: Int32[1, 1, 3, 3] primitive_lattice: [3.111 -1.5555 0.0; 0.0 2.6942050311733885 0.0; 0.0 0.0 4.988] mapping_to_primitive: Int32[1, 2, 3, 4] n_std_atoms: 4 std_lattice: [3.111 -1.5555 0.0; 0.0 2.6942050311733885 0.0; 0.0 0.0 4.988] std_types: Int32[1, 1, 2, 2] std_positions: StaticArraysCore.SVector{3, Float64}[[0.3333333333333333, 0.6666666666666666, 0.0], [0.6666666666666667, 0.3333333333333333, 0.5], [0.3333333333333333, 0.6666666666666666, 0.6181], [0.6666666666666667, 0.3333333333333333, 0.1181000000000001]] std_rotation_matrix: [1.0 0.0 0.0; 0.0 1.0 0.0; 0.0 0.0 1.0] std_mapping_to_primitive: Int32[1, 2, 3, 4] pointgroup_symbol: 6mm
julia> dataset.primitive_lattice3×3 Lattice{Float64} 3.111 -1.5555 0.0 0.0 2.6942050311733885 0.0 0.0 0.0 4.988
julia> dataset.std_lattice3×3 Lattice{Float64} 3.111 -1.5555 0.0 0.0 2.6942050311733885 0.0 0.0 0.0 4.988
julia> std_lattice_before_idealization = convert(Matrix{Float64}, Lattice(cell)) * inv(dataset.transformation_matrix)3×3 Matrix{Float64}: 3.111 -1.5555 0.0 0.0 2.69421 0.0 0.0 0.0 4.988