matador.crystal package

This module contains functionality for wrapping matador data dictionaries as Crystal objects, which possess several convenience methods.

class matador.crystal.Crystal(doc, voronoi=False, network_kwargs=None, mutable=False)[source]

Bases: DataContainer

Class that wraps the MongoDB document, providing useful interfaces for cell manipulation and validation.

elems

list of present elements in the crystal, sorted in alphabetical order.

Type:

list of str

cell

the unit cell constructed from lattice vectors.

Type:

UnitCell

Initialise Crystal object from matador document with Site list and any additional abstractions, e.g. voronoi or CrystalGraph.

Parameters:

doc (dict) – document containing structural information, minimal requirement is for atom_types, one of lattice_abc or lattice_cart, and one of positions_frac or positions_abs to be present.

Keyword Arguments:
  • voronoi (bool) – whether to compute Voronoi substructure for each site

  • network_kwargs (dict) – keywords to pass to the CrystalGraph initialiser

update(data)[source]

Update the underlying self._data dictionary with the passed data.

print_sites()[source]
set_positions(new_positions, fractional=True)[source]
property atom_types: List[str]

Return list of atom types.

property num_atoms: int

Return number of atoms in structure.

property num_elements: int

Return number of species in the structure.

property positions_frac: List[List[float]]

Return list of fractional positions.

property positions_abs: List[List[float]]

Return list of absolute Cartesian positions.

property site_occupancies: List[float]

Return a list of site occupancies.

property stoichiometry: List[List[Union[str, float]]]

a list of two-member lists containing element symbol and number of atoms per formula unit, sorted in alphabetical order by element symbol).

Type:

Return stoichiometry in matador format

property num_fu: int
property concentration: List[float]

Return concentration of each species in stoichiometry.

get_concentration(elements=None, include_end=True) List[float][source]

Return concentration of each species in stoichiometry, and use this value for future calls to crystal.concentration.

property formula: str

Returns chemical formula of structure.

property formula_tex: str

Returns chemical formula of structure in LaTeX format.

property formula_unicode: str

Returns a unicode string for the chemical formula.

property cell_volume: float

Returns cell volume in ų.

property lattice_cart: Tuple[Tuple[float, float, float], Tuple[float, float, float], Tuple[float, float, float]]

The Cartesian lattice vectors as a tuple.

property lattice_abc: Tuple[Tuple[float, float, float], Tuple[float, float, float]]

Lattice parameters as a tuple.

property bond_lengths: List[Tuple[Tuple[str, str], float]]

Returns a list of ((species_A, species_B), bond_length)), sorted by bond length, computed from the network structure of the crystal (i.e. first coordination sphere).

property voronoi_substructure

Returns, or calculates if not present, the Voronoi substructure of crystal.

property space_group: str

Return the space group symbol at the last-used symprec.

property space_group_tex: str

Returns the space group symbol at the last-used symprec, formatted for LaTeX.

The HM symbol is surrounded by $ and with ‘-’ replaced with overbars for rotoinversion axes, e.g., Fm-3m -> $Fmbar{3}m$.

get_space_group(symprec=0.01) str[source]

Return the space group of the structure at the desired symprec. Stores the space group in a dictionary self._space_group under symprec keys. Updates self._data['space_group'] and self._data['symprec'] with the last value calculated.

Keyword Arguments:

symprec (float) – spglib symmetry tolerance.

property pdf

Returns a PDF object (pair distribution function) for the structure, calculated with default PDF settings.

calculate_pdf(**kwargs)[source]

Calculate and set the PDF with the passed parameters.

property pxrd

Returns a PXRD object (powder xray diffraction) containing the XRD pattern for the structure.

calculate_pxrd(**kwargs)[source]

Compute and set PXRD with the passed parameters.

property coordination_stats

Returns stastistics on the coordination of each site, as computed from Voronoi decomposition.

property ase_atoms

Returns an ASE Atoms representation of the crystal.

property pmg_structure

Returns the pymatgen structure representation of the crystal.

classmethod from_ase(atoms)[source]
classmethod from_pmg(pmg)[source]
property coordination_lists

Returns a dictionary containing pairs of elements and the list of coordination numbers per site.

property unique_sites

Return unique sites using Voronoi decomposition.

property network

Returns/constructs a CrystalGraph object of the structure.

property network_stats

Return network-calculated coordnation stats.

property bonding_stats

Return network-calculated bonding stats.

Returns:

sorted dictionary with root atom as keys and bond

information as values.

Return type:

dict

draw_network(layout=None)[source]

Draw the CrystalGraph network.

supercell(extension=None, target=None)[source]

Returns a supercell of this crystal with the specified extension.

standardized(primitive=False)[source]

Returns a standardized cell for this crystal, optionally the primitive cell.

Submodules

matador.crystal.crystal module

This submodule implements the Crystal class, a wrapper to the raw dictionary stored in MongoDB that allows for validation, manipulation and analysis of the lattice.

class matador.crystal.crystal.UnitCell(lattice)[source]

Bases: object

This class describes a 3D periodic unit cell by its Cartesian lattice vectors or lattice parameters, in Å.

Initialise the cell from either Cartesian lattice vectors [[x1, y1, z1], [x2, y2, z2], [x3, y3, z3]], or lattice parameters [[a, b, c], [alpha, beta, gamma]].

Parameters:

lattice (list or numpy.ndarray) – either Cartesian lattice vectors or lattice parameters, stored as lists or a numpy arrays.

property lattice_cart: Tuple[Tuple[float, float, float], Tuple[float, float, float], Tuple[float, float, float]]

The Cartesian lattice vectors as a tuple.

property lattice_abc: Tuple[Tuple[float, float, float], Tuple[float, float, float]]

Lattice parameters as a tuple.

property lengths: Tuple[float, float, float]

Lattice vector lengths.

property recip_lattice_cart: List[List[float]]
property angles: Tuple[float, float, float]

Lattice vector angles.

property volume: float

The cell volume in ų.

class matador.crystal.crystal.Crystal(doc, voronoi=False, network_kwargs=None, mutable=False)[source]

Bases: DataContainer

Class that wraps the MongoDB document, providing useful interfaces for cell manipulation and validation.

elems

list of present elements in the crystal, sorted in alphabetical order.

Type:

list of str

cell

the unit cell constructed from lattice vectors.

Type:

UnitCell

Initialise Crystal object from matador document with Site list and any additional abstractions, e.g. voronoi or CrystalGraph.

Parameters:

doc (dict) – document containing structural information, minimal requirement is for atom_types, one of lattice_abc or lattice_cart, and one of positions_frac or positions_abs to be present.

Keyword Arguments:
  • voronoi (bool) – whether to compute Voronoi substructure for each site

  • network_kwargs (dict) – keywords to pass to the CrystalGraph initialiser

update(data)[source]

Update the underlying self._data dictionary with the passed data.

print_sites()[source]
set_positions(new_positions, fractional=True)[source]
property atom_types: List[str]

Return list of atom types.

property num_atoms: int

Return number of atoms in structure.

property num_elements: int

Return number of species in the structure.

property positions_frac: List[List[float]]

Return list of fractional positions.

property positions_abs: List[List[float]]

Return list of absolute Cartesian positions.

property site_occupancies: List[float]

Return a list of site occupancies.

property stoichiometry: List[List[Union[str, float]]]

a list of two-member lists containing element symbol and number of atoms per formula unit, sorted in alphabetical order by element symbol).

Type:

Return stoichiometry in matador format

property num_fu: int
property concentration: List[float]

Return concentration of each species in stoichiometry.

get_concentration(elements=None, include_end=True) List[float][source]

Return concentration of each species in stoichiometry, and use this value for future calls to crystal.concentration.

property formula: str

Returns chemical formula of structure.

property formula_tex: str

Returns chemical formula of structure in LaTeX format.

property formula_unicode: str

Returns a unicode string for the chemical formula.

property cell_volume: float

Returns cell volume in ų.

property lattice_cart: Tuple[Tuple[float, float, float], Tuple[float, float, float], Tuple[float, float, float]]

The Cartesian lattice vectors as a tuple.

property lattice_abc: Tuple[Tuple[float, float, float], Tuple[float, float, float]]

Lattice parameters as a tuple.

property bond_lengths: List[Tuple[Tuple[str, str], float]]

Returns a list of ((species_A, species_B), bond_length)), sorted by bond length, computed from the network structure of the crystal (i.e. first coordination sphere).

property voronoi_substructure

Returns, or calculates if not present, the Voronoi substructure of crystal.

property space_group: str

Return the space group symbol at the last-used symprec.

property space_group_tex: str

Returns the space group symbol at the last-used symprec, formatted for LaTeX.

The HM symbol is surrounded by $ and with ‘-’ replaced with overbars for rotoinversion axes, e.g., Fm-3m -> $Fmbar{3}m$.

get_space_group(symprec=0.01) str[source]

Return the space group of the structure at the desired symprec. Stores the space group in a dictionary self._space_group under symprec keys. Updates self._data['space_group'] and self._data['symprec'] with the last value calculated.

Keyword Arguments:

symprec (float) – spglib symmetry tolerance.

property pdf

Returns a PDF object (pair distribution function) for the structure, calculated with default PDF settings.

calculate_pdf(**kwargs)[source]

Calculate and set the PDF with the passed parameters.

property pxrd

Returns a PXRD object (powder xray diffraction) containing the XRD pattern for the structure.

calculate_pxrd(**kwargs)[source]

Compute and set PXRD with the passed parameters.

property coordination_stats

Returns stastistics on the coordination of each site, as computed from Voronoi decomposition.

property ase_atoms

Returns an ASE Atoms representation of the crystal.

property pmg_structure

Returns the pymatgen structure representation of the crystal.

classmethod from_ase(atoms)[source]
classmethod from_pmg(pmg)[source]
property coordination_lists

Returns a dictionary containing pairs of elements and the list of coordination numbers per site.

property unique_sites

Return unique sites using Voronoi decomposition.

property network

Returns/constructs a CrystalGraph object of the structure.

property network_stats

Return network-calculated coordnation stats.

property bonding_stats

Return network-calculated bonding stats.

Returns:

sorted dictionary with root atom as keys and bond

information as values.

Return type:

dict

draw_network(layout=None)[source]

Draw the CrystalGraph network.

supercell(extension=None, target=None)[source]

Returns a supercell of this crystal with the specified extension.

standardized(primitive=False)[source]

Returns a standardized cell for this crystal, optionally the primitive cell.

matador.crystal.crystal_site module

This submodule implements the Site class for handling atomic sites.

class matador.crystal.crystal_site.Site(species: str, position: list, lattice, position_unit='fractional', mutable: bool = False, **site_data)[source]

Bases: DataContainer

The Site class contains a description of an individual site within a 3D periodic Crystal.

Initialise a Site object from its species, position and a reference to the lattice it exists in. Any other keys will be made available as site-level values.

get(key, default=None)[source]

Overload dictionary.get() method.

Parameters:

key (str) – key to try and obtain.

Keyword Arguments:

default – return value raise if key is not present.

Returns:

_data[key] if it exists, else None.

set_position(position, units)[source]
property coords
property coords_cartesian
property species
property lattice
property occupancy
property coordination
displacement_between_sites(other_site)[source]
distance_between_sites(other_site)[source]

matador.crystal.elastic module

This submodule contains functionality to fit E(V) equations of state to ab initio data, primarily to calculate bulk moduli.

matador.crystal.elastic.get_equation_of_state(seed, plot=False)[source]

Extract E(V) data from CASTEP files and perform fits for the equation of state and bulk moduli.

Parameters:

seed (str/dict) – filename or scraped dictionary to fit.

Keyword Arguments:

plot (bool) – plot the fitted EoS.

class matador.crystal.elastic.EquationOfState(E_0, V_0)[source]

Bases: object

Abstract class for E(V) isothermal equations of state. Used to perform least squares fitting to arbitrary functional form.

E_0

equilbirium lattice energy.

Type:

float

V_0

equilibrium lattice volume.

Type:

float

popt

fitting parameters for particular EOS.

Type:

list of float

p0

initial guess at fitting parameters.

Type:

list of float

rrmsd

relative root mean square deviation of fit, as defined by [2].

Type:

float

Set up EOS ready for fit.

Parameters:
  • E_0 (float) – equilibrium energy.

  • V_0 (float) – equilibrium volume.

fit(volumes, energies)[source]

Perform a least squares fit on the volumes and energies.

Parameters:
  • volumes (list of float) – calculated volumes.

  • energies (list of float) – calculated energies.

evaluate(V, B, C)[source]

Evaluate the EOS.

Parameters:
  • V (list of float) – volumes at which to test.

  • B (float) – the first fitting parameter defined in [2].

  • C (float) – the second fitting parameter defined in [2].

residual(guess, V, E)[source]

Calculate the resdiual of the current fit.

Parameters:
  • guess (list of float) – interim fitting parameters.

  • V (list of float) – volumes at which to test.

  • E (list of float) – energies at the volumes V.

property bulk_modulus

Returns the bulk modulus predicted by the fit.

property bulk_modulus_err

Returns the estimated error on the bulk modulus.

property fit_parameters

Return list of final fitting parameters.

class matador.crystal.elastic.BirchMurnaghanEulerianEOS(E_0, V_0)[source]

Bases: EquationOfState

Implements the 3rd order Birch-Murnaghan EOS [1] for given data, as provided by Ref. [2] in the Eulerian frame.

[1]. Francis Birch, Phys. Rev. 71 809 (1947)

DOI: 10.1103/PhysRev.71.809.

[2]. K. Latimer, S. Dwaraknath, K. Mathew, D. Winston, K. A. Persson,

npj Comput. Mater. 2018 41 2018, 4, 40. DOI: 10.1038/s41524-018-0091-x.

Set up EOS ready for fit.

Parameters:
  • E_0 (float) – equilibrium energy.

  • V_0 (float) – equilibrium volume.

evaluate(V, B, C)[source]

Evaluate the EOS.

Parameters:
  • V (list of float) – volumes at which to test.

  • B (float) – the first fitting parameter defined in [2].

  • C (float) – the second fitting parameter defined in [2].

get_bulk_modulus()[source]

Return the bulk modulus of this particular fit.

class matador.crystal.elastic.PoirerTarantolaEOS(E_0, V_0)[source]

Bases: EquationOfState

Implements the logarithmic Poirer-Tarantola [3] for given data, as provided by Ref. [2].

[2]. K. Latimer, S. Dwaraknath, K. Mathew, D. Winston, K. A. Persson,

npj Comput. Mater. 2018 41 2018, 4, 40. DOI: 10.1038/s41524-018-0091-x.

[3]. J. Poirer, A. Tarantola, Phys. Earth Planet. Inter. 109 1-8 (1998).

Set up EOS ready for fit.

Parameters:
  • E_0 (float) – equilibrium energy.

  • V_0 (float) – equilibrium volume.

evaluate(V, B, C)[source]

Evaluate the EOS.

Parameters:
  • V (list of float) – volumes at which to test.

  • B (float) – the first fitting parameter defined in [2].

  • C (float) – the second fitting parameter defined in [2].

get_bulk_modulus()[source]

Return the bulk modulus of this particular fit.

class matador.crystal.elastic.TaitEOS(E_0, V_0)[source]

Bases: EquationOfState

Implements the exponential Tait EOS [4] for given data, as provided by Ref. [2].

[2]. K. Latimer, S. Dwaraknath, K. Mathew, D. Winston, K. A. Persson,

npj Comput. Mater. 2018 41 2018, 4, 40. DOI: 10.1038/s41524-018-0091-x.

[4]. J. H. Dymond, R. Malhotra, Int. J. Thermophys. 9 941-951 (1988).

Set up EOS ready for fit.

Parameters:
  • E_0 (float) – equilibrium energy.

  • V_0 (float) – equilibrium volume.

evaluate(V, B, C)[source]

Evaluate the EOS.

Parameters:
  • V (list of float) – volumes at which to test.

  • B (float) – the first fitting parameter defined in [2].

  • C (float) – the second fitting parameter defined in [2].

get_bulk_modulus()[source]

Return the bulk modulus of this particular fit.

matador.crystal.elastic.plot_volume_curve(probes_all, labels, curves, volumes, energies, seed)[source]

Plot all fitted EOS curves for this structure and save as a PDF.

matador.crystal.network module

This file implements turning matador Crystal objects into CrystalGraph objects.

class matador.crystal.network.CrystalGraph(structure=None, graph=None, coordination_cutoff=1.1, bond_tolerance=1e+20, num_images=1, debug=False, separate_images=False, delete_one_way_bonds=False, max_bond_length=5)[source]

Bases: MultiDiGraph

Create networkx.MultiDiGraph object with extra functionality for atomic networks.

Keyword Arguments:
  • structure (matador.Crystal) – crystal structure to network-ify

  • graph (nx.MultiDiGraph) – initialise from graph

  • coordination_cutoff (float) – max multiplier of first coordination sphere for edge drawing

  • num_images (int) – number of periodic images to include in each direction

  • separate_images (bool) – whether or not to include image atoms as new nodes

get_strongly_connected_component_subgraphs(delete_one_way_bonds=True)[source]

Return generator of strongly-connected subgraphs in CrystalGraph format.

get_communities(graph=None, **louvain_kwargs)[source]

Return list of community subgraphs in CrystalGraph format.

remove_directionality(graph=None)[source]
set_unique_subgraphs(method='community')[source]

Filter strongly connected component subgraphs for isomorphism with others inside CrystalGraph. Sets self.unique_subgraph to a set of such subgraphs.

get_bonds_per_atom()[source]
matador.crystal.network.node_match(n1, n2)[source]
matador.crystal.network.get_unique_subgraphs(subgraphs)[source]

Filter strongly connected component subgraphs for isomorphism with others.

Input:

subgraphs: list(CrystalGraph), list of subgraph objects to filter
Returns:

set(CrystalGraph), set of unique subgraphs

Return type:

unique_subgraphs

matador.crystal.network.are_graphs_the_same(g1, g2, edge_match=None)[source]
matador.crystal.network.draw_network(structure, layout=None, edge_labels=False, node_index=False, curved_edges=True, node_colour='elem', partition=None, ax=None)[source]