Structure Module

Contains a minimalist SmactStructure class for simple crystal structure and chemical composition representation and manipulation.

Minimalist structure representation for comprehensible manipulation.

class smact.structure_prediction.structure.SmactStructure[source]

Bases: object

SMACT implementation inspired by pymatgen Structure class.

Handles basic structural and compositional information for a compound. Includes a lossless POSCAR-style specification for storing structures, allowing structures to be stored in files or databases, or to be pulled from the Materials Project.

species

A list of tuples describing the composition of the structure, stored as (element, oxidation, stoichiometry). The list is sorted alphabetically based on element symbol, and identical elements are sorted with highest charge first.

lattice_mat

A numpy 3x3 array containing the lattice vectors.

sites

A dictionary of {species: coords}, where species is a string representation of the species and coords is a list of position vectors, given as lists of length 3. For example:

>>> s = SmactStructure.from_file('tests/files/NaCl.txt')
>>> s.sites
{'Cl1-': [[2.323624165, 1.643050405, 4.02463512]], 'Na1+': [[0.0, 0.0, 0.0]]}
lattice_param

The lattice parameter.

as_poscar() str[source]

Represent the structure as a POSCAR file compatible with VASP5.

The POSCAR format adopted is as follows:

The first line contains the species’ names separated by a whitespace. The second through fourth line, inclusive, contain the lattice matrix: each line contains a lattice vector, with elements separated by a whitespace. The fifth line contains the elements’ names separated by a whitespace. If more than one oxidation state exists for an element, the element appears multiple times; once for each oxidation state. The sixth line is the string ‘Cartesian’. The seventh line onwards are the Cartesian coordinates of each site, separated by a whitespace. In addition, at the end of each line is the species’ name, separated by a whitespace.

For examples of this format, see the text files under tests/files.

Returns

POSCAR-style representation of the structure.

Return type

str

composition() str[source]

Generate a key that describes the composition.

Key format is ‘{element}_{stoichiometry}_{charge}{sign}’ with no delimiter, sans brackets. Species are ordered as stored within the structure, see SmactStructure.

Returns

Key describing constituent species.

Examples

>>> s = SmactStructure.from_file('tests/files/CaTiO3.txt')
>>> print(s.composition())
Ca_1_2+O_3_2-Ti_1_4+
static from_file(fname: str)[source]

Create SmactStructure from a POSCAR file.

Parameters

fname – The name of the POSCAR file. See as_poscar() for format specification.

Returns

SmactStructure

static from_mp(species: List[Union[Tuple[str, int, int], Tuple[Species, int]]], api_key: str, determine_oxi: str = 'BV')[source]

Create a SmactStructure using the first Materials Project entry for a composition.

Parameters
  • species – See __init__().

  • determine_oxi (str) – The method to determine the assignments oxidation states in the structure. Options are ‘BV’, ‘comp_ICSD’,’both’ for determining the oxidation states by bond valence, ICSD statistics or trial both sequentially, respectively.

  • api_key – A www.materialsproject.org API key.

Returns

SmactStructure

static from_poscar(poscar: str)[source]

Create SmactStructure from a POSCAR string.

Parameters

poscar – A SMACT-formatted POSCAR string. See as_poscar() for format specification.

Returns

SmactStructure

static from_py_struct()[source]

Create a SmactStructure from a pymatgen Structure object.

Parameters
  • structure – A pymatgen Structure.

  • determine_oxi (str) – The method to determine the assignments oxidation states in the structure. Options are ‘BV’, ‘comp_ICSD’,’both’ for determining the oxidation states by bond valence, ICSD statistics or trial both sequentially, respectively.

Returns

SmactStructure

get_spec_strs() List[str][source]

Get string representations of the constituent species.

Returns

A list of strings, formatted as ‘{element}{charge}{sign}’.

Examples

>>> s = SmactStructure.from_file('tests/files/CaTiO3.txt')
>>> s.get_spec_strs()
['Ca2+', 'O2-', 'Ti4+']
has_species(species: Tuple[str, int]) bool[source]

Determine whether a given species is in the structure.