CellBase.jl

Documentation for CellBase.jl

What is package does

Provide a basic low level interface for storing and handling crystal structure data. The target application is for small and periodic cells, with the emphasis on both ease to use as well as performance.

Guides

Lattice type

The Lattice type represents the lattice. The lattice vectors are stored as a matrix of column vectors. For example, the first lattice vector should be accessed as cellmat(lattice)[:, 1].

Updating an existing `Lattice`

The reciprocal lattice is also stored for quick access, however, this means that whenever the matrix is directly modified, the store reciprocal lattice vectors should be updated as well. Hence, the set_cellmat! function should be used for updating the cell matrix, which does this is automatically.

CellBase.LatticeType

Lattice{T}

The type represents a lattice in 3D space formed by three column lattice vectors.

CellBase.LatticeMethod
Lattice(matrix::Matrix{T}) where T

Construct a Lattice from a matrix of column vectors.

CellBase.LatticeMethod
Lattice(a::T, b::T, c::T, α::T, β::T, γ::T) where T

Construct a Lattice from lattice parameters.

CellBase.LatticeMethod
Lattice(a::T, b::T, c::T) where T <: Real

Construct a Lattice with orthogonal lattice vectors.

CellBase.LatticeMethod
Lattice(va::Vector{T}, vb::Vector{T}, vc::Vector{T}) where T

Construct a Lattice from three lattice vectors.

CellBase.LatticeMethod
Lattice(cellpar::Vector{T}) where T

Construct a Lattice from a six-vector of lattice parameters.

CellBase.cellmatMethod
cellmat(lattice::Lattice)

Return the matrix of column vectors.

CellBase.cellparMethod
cellpar(lattice::Lattice)

Return the lattice parameters as a six-vector.

CellBase.micMethod
mic(l::Lattice, vec::AbstractVecOrMat)

Compute the minimum-image convention representation of a series of displacement vectors.

CellBase.mic_naiveMethod

Compute naive MIC representation of the vector(s)

Not safe for skewed cells. Requires

norm(length) < 0.5 * min(cellpar(l)[1:3])

See also:

  • W. Smith, "The Minimum Image Convention in Non-Cubic MD Cells", 1989, http://citeseerx.ist.psu.edu/viewdoc/summary?doi=10.1.1.57.1696.
CellBase.mic_naiveMethod
 mic_safe(l::Lattice, v::AbstractVector)

Compute MIC representation of vectors using a safe approach based on Minkowski reduced cell

CellBase.mic_safeMethod
 mic_safe(l::Lattice, v::AbstractMatrix)

Compute MIC representation of vectors using a safe approach based on Minkowski reduced cell

CellBase.mic_safeMethod
 mic_safe(l::Lattice, v::AbstractVector)

Compute MIC representation of vectors using a safe approach based on Minkowski reduced cell

CellBase.random_vec_in_cellMethod
random_vec_in_cell(cell::Matrix{T}; scales::Vector=ones(size(cell)[1])) where T

Get an random vector within a unit cell. The cell is a matrix made of column vectors.

CellBase.rec_cellmatMethod
rec_cellmat(l::Lattice)

Returns the matrix of reciprocal lattice vectors (without the $2\pi$ factor).

CellBase.reciprocalMethod
reciprocal(l::Lattice)

Returns the matrix of reciprocal lattice vectors (without the $2\pi$ factor).

CellBase.set_cellmat!Method
set_cellmat!(lattice::Lattice, mat)

Set the cell matrix and update the reciprocal cell matrix.

CellBase.update_rec!Method
update_rec!(l::Lattice)

Update the reciprocal matrix - this should be called everytime cell is changed.

CellBase.wrap_positionsMethod
wrap_positions(l::Lattice, v::AbstractVecOrMat)

Wrap positions back into the box defined by the cell vectors.

CellBase.max_shiftsMethod

Compute the require periodicities required to fill a cut off sphere of a certain radius

CellBase.shift_indicesMethod
shift_indices(lattice::AbstractArray, rc::Real)

Compute a vector containing shift indices

CellBase.shift_vectorsFunction
shift_vectors(lattice::AbstractMatrix, shift1, shift2, shift3)

Compute the shift vectors needed to include all parts of the lattice within a cut off radius.

Cell type

The Cell type represent a crystal structure combining the positions of atoms and the lattice vectors. It basically combines the Lattice, the atomic positions and their identities. Internally, the formere is stored as an matrix of column vectors of the absolute positions (cartesian space).

For flexibility, any additional array like data can be stored in the arrays property which is a dictionary with Symbol keys. Additional metadata can be also be stored under metadata, which is a Dict{Any, Any}.

CellBase.CellType

A Cell represents a periodic structure in three-dimensional space.

Defined as:

mutable struct Cell{T}
    lattice::Lattice{T}                 # Lattice of the structure
    symbols::Vector{Symbol}
    positions::Matrix{T}
    arrays::Dict{Symbol, Any}        # Any additional arrays
    metadata::Dict{Symbol, Any}
end
CellBase.CellMethod
Cell(lat::Lattice, numbers::Vector{Int}, positions::Vector)

Constructure the Cell type from lattice, positions and numbers

CellBase.CellMethod
Cell(l::Lattice, symbols, positions) where T

Construct a Cell type from arrays

CellBase.CellMethod
Cell(lat::Lattice, numbers::Vector{Int}, positions)

Constructure the Cell type from lattice, positions and numbers

Base.sort!Method
Base.sort!(cell::Cell)

In-place sort a Cell with the positions sorted by the species kinds.

Note

This will in-place-update all underlying Array - use with caution where the data is shared between multiple instances.

Base.sortMethod
Base.sort(cell::Cell)

Sort a Cell with the positions sorted by the species kinds.

CellBase._distance_matrix_micMethod
_distance_matrix_mic(cell::Cell)

Compute the distance matrix for the given structure, using the minimum image convention (or not).

Note the returned matrix does not expand the cell. The distance matrix cannot be safety used for obtaining the minimum separations. For example, a structure with a single atom would be a distance matrix containing only zero.

CellBase._distance_matrix_no_micMethod
_distance_matrix_no_mic(structure::Cell)

Compute the distnace matrix without minimum image convention, e.g. the PBC is not respected.

CellBase.arrayMethod
array(structure::Cell, arrayname::Symbol)

Return the additional array stored in the Cell object.

CellBase.arraynamesMethod
arraynames(structure::Cell)

Return the names of additional arrays.

CellBase.attachmetadata!Method
attachmetadata!(structure::Cell, metadata::Dict)

Replace metadata with an existing dictionary.

CellBase.cellmatMethod
cellmat(structure::Cell)

Return the matrix of lattice vectors.

CellBase.cellparMethod
cellpar(structure::Cell)

Return the lattice parameters.

CellBase.check_minsepMethod
check_minsep(structure::Cell, minsep::Dict)

Check if the minimum separation constraints are satisfied. Minimum separations are supplied as an dictionary where the global version under the :global key. To supply the minimum separations between A and B, pass Dict((:A=>:B)=>1.0). Return true or false.

CellBase.clipMethod
clip(s::Cell, mask::AbstractVector)

Clip a structure with a given indexing array

CellBase.compute_minsep_matMethod
minsep_matrix(structure::Cell, minsep::Dict)

Initialise the minimum separation matrix and species mapping. Returns the unique species, integer indexed species and the minimum separation matrix.

CellBase.distance_squared_betweenMethod
distance_squared_between(posmat::Matrix, i, j, svec::Matrix, ishift)

Return the squared distance between two positions stored in a matrix and shift vector.

CellBase.fingerprintMethod
fingerprint(s::Cell; dmat=distance_matrix(s), weighted=true, cut_bl=3.0)

Computed the fingerprint vector based on simple sorted pair-wise distances. NOTE: Does not work for single atom cell!!

CellBase.fingerprint_distanceMethod
fingerprint_distance(f1::AbstractVector, f2::AbstractVector;lim=Inf)

Compute the deviation between two finger print vectors

Comparison is truncated by the size of the shortest vector of the two, or by the lim key word.

CellBase.get_cellmatMethod
get_cellmat(structure::Cell)

Return the matrix of lattice vectors(copy).

CellBase.get_positionsMethod
get_positions(cell::Cell)

Return a copy of the positions (cartesian coordinates) of the atoms in a structure.

CellBase.make_supercellMethod
make_supercell(structure::Cell, a, b, c)

Make a supercell

Currently only work with diagonal transform matrices. TODO: Write function for the general cases....

CellBase.metadataMethod
metadata(structure::Cell)

Return the metadata dictionary.

CellBase.natomsFunction
natoms(cell::Cell)

Return number of atoms in a structure.

CellBase.nionsMethod
nions(cell::Cell)

Return number of atoms in a structure.

CellBase.num_fuMethod
num_fu(structure::Cell)

Return the number of formula units.

CellBase.positionsMethod
positions(cell::Cell)

Return positions (cartesian coordinates) of the atoms in a structure.

CellBase.rattle!Method
rattle!(cell::Cell, amp)

Rattle the positions of the cell for a given maximum amplitude (uniform distribution).

CellBase.set_cellmat!Method
set_cellmat!(cell::Cell, mat;scale_positions=true)

Update the Lattice with a new matrix of lattice vectors. Scale of the existing postions if needed.

CellBase.speciesMethod
species(structure::Cell)

Return a Vector of species names.

CellBase.specindexMethod
specindex(structure::Cell)

Return the unique species and integer based indices for each atom.

CellBase.sposarrayMethod
sposarray(cell::Cell)

Return the positions as a Vector of static arrays. The returned array can provide improved performance for certain type of operations.

CellBase.volumeMethod
volume(structure::Cell)

Return the volume of the cell.

CellBase.wrap!Method
wrap!(vec::AbstractVector, l::Lattice)

Wrap a vector back to the periodic box defined by the lattice vectors.

CellBase.wrap!Method
wrap!(cell::Cell)

Wrap an atom outside of the lattice back into the box defined by the lattice vectors.

Spglib.jl interface

The routine in Spglib.jl (which wraps the spglib libraray) can be used for finding symmetry and performing reduction of the Cell type.

CellBase.CellMethod
Cell(cell::SCell)

Return a Cell object from Spglib.Cell where the SCell.types have been changed to 1-based index. This is useful when the Spglib._expand_cell is called where the types returned is re-indexed.

CellBase.CellMethod
Cell(cell::SCell)

Return a Cell object from Spglib.Cell.

CellBase.SCellMethod
SCell(cell::Cell)

Construct Spglib.Cell from Cell type.

CellBase.niggli_reduce_cellFunction
niggli_reduce_cell(cell::Cell; wrap_pos=true)

Apply niggli reduction to the lattice using Spglib

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

Apply niggli reduction to the lattice using Spglib. The positions are not wrapped.

CellBase.@extend_scellMacro

Macro for extending the Spglib methods.

Usage:

@extend_scell get_dataset

will allow the get_dataset method of Spglib to be used for Cell type.

CellBase.@extend_scell_roundtripMacro

Macro for extending the Spglib methods and convert returned Spglib.Cell to Cell.

Usage:

@extend_scell_roundtrip standardize_cell

will allow the standardize_cell method to be used and the returned Spglib.Cell is converted to Cell. This assumes the only changes made is on the Lattice and there is no change in the atomic positions/orders.

CellBase.@extend_scell_roundtrip_newMacro

Macro for extending the Spglib methods and convert returned Spglib.Cell to Cell.

Usage:

@extend_scell_roundtrip_new find_primitive

will allow the standardize_cell method to be used and the returned Spglib.Cell is converted to Cell.

Misc utils

Miscellaneous utility functions.

CellBase.vec2cellparMethod

Convert cell vectors to cell parameters

Returns an static array of the cell parameters

Neighbour lists

Type and functions for handling neighbour lists.

CellBase.ExtendedPointArrayMethod
ExtendedPointArray(cell::Cell, rcut)

Constructed an ExtendedPointArray from a given structure. Implicitly, the positions are wrapped inside the unit cell, even if the actual in the original Cell is outside the unit cell. This ensures the correct neighbour list begin constructed.

CellBase.NeighbourListMethod
NeighbourList(ea::ExtendedPointArray, rcut, nmax=1000; savevec=false)

Construct a NeighbourList from an extended point array for the points in the original cell

CellBase._need_rebuildMethod
_need_rebuild(nl::NeighbourList)

Check if a full rebuild is needed for the NeighbourList if the skin is used and if any atom has move more than the skin.

CellBase.eachneighbourMethod

Iterate the neighbours of a site in the original cell. Returns a tuple of (originalindex, extendedindex, distance) for each iteration

CellBase.eachneighbourvectorMethod

Iterate the neighbours of a site in the original cell. Returns a tuple of (originalindex, extendedindex, distance, vector) for each iteration

CellBase.get_neighbourMethod
get_neighbour(nl::NeighbourList, iorig::Int, idx::Int)

Return a single neighbour.

CellBase.get_neighboursMethod
get_neighbours(nl::NeighbourList, iorig::Int)

Return a Vector of Neighbour objects

CellBase.rebuild!Method
rebuild!(ea::ExtendedPointArray, cell)

Rebuild the ExtendedPointArray for an existing cell

CellBase.rebuild!Method
rebuild(nl::NeighbourList, cell::Cell)

Perform a full rebuild of the NeighbourList with the latest geometry of the cell

CellBase.rebuild!Method
rebuild!(nl::NeighbourList, ea::ExtendedPointArray)

Perform a full rebuild of the neighbour list from scratch for a given ExtendedPointArray. Extended the neighbour storage space if necessary.

CellBase.update!Method
update!(nl::NeighbourList, cell::Cell)

Update the NeighbourList with the latest geometry of the Cell. No rebuilding is performed.

CellBase.update!Method

Update the calculated distances/vectors but do not rebuild the whole neighbour list

CellBase.update_nmax!Method
increase_nmax!(nl::NeighbourList, nmax)

Increase the maximum number of neighbours storable in the neighbour list

IO

Code for read/writing file of crystal structures.

CellBase.write_resFunction
write_res(fname::AbstractString, structure::Cell, mode="w")

Write out SHELX format data to a file

CellBase.write_resMethod
write_res(io::IO, structure::Cell)

Write out SHELX format data including additional information stored in structure.metadata.

Supported keys: :label, :pressure, :volume, :enthalpy, :spin, :abs_spin, :symm, :flag1, flag2, flag3.

The following fields are also treated as the REM line:

  • info: If a Dict is supplied the key-value pairs will be written.
  • comments: Any additional comments to be written. Must be supplied as a Vector{<:AbstractString}.

The composition of the structure will be written as REM Composition: <E1> <N1> <E2> <N2>.

CellBase.read_xyzMethod

Read XYZ file

NOTE: This is a pure julia implementation which may not work with all cases. In the future, better to use ExtXYZ.jl interface.

CellBase.SheapIO.run_sheapMethod
run_sheap(vecs, opt::SheapOptions;metadata=repeat([SheapMetadata()], length(vecs)), show_stderr=true)

Run SHEAP for a iterator of vectors with the given options.