pymatgen.analysis.local_env module

class JMolNN(tol=0.001, el_radius_updates=None)[source]

Bases: pymatgen.analysis.local_env.NearNeighbors

Determine near-neighbor sites and coordination number using an emulation of JMol’s default autoBond() algorithm. This version of the algorithm does not take into account any information regarding known charge states.

Parameters:
  • tol (float) – tolerance parameter for bond determination (default: 1E-3).
  • el_radius_updates – (dict) symbol->float to override default atomic radii table values
get_max_bond_distance(el1_sym, el2_sym, constant=0.56)[source]

Use JMol algorithm to determine bond length from atomic parameters :param el1_sym: (str) symbol of atom 1 :param el2_sym: (str) symbol of atom 2 :param constant: (float) factor to tune model

Returns: (float) max bond length

get_nn_info(structure, n)[source]

Get all near-neighbor sites as well as the associated image locations and weights of the site with index n using the bond identification algorithm underlying JMol.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine near neighbors.
Returns:

tuples, each one

of which represents a neighbor site, its image location, and its weight.

Return type:

siw (list of tuples (Site, array, float))

class MinimumDistanceNN(tol=0.1, cutoff=10.0)[source]

Bases: pymatgen.analysis.local_env.NearNeighbors

Determine near-neighbor sites and coordination number using the nearest neighbor(s) at distance, d_min, plus all neighbors within a distance (1 + delta) * d_min, where delta is a (relative) distance tolerance parameter.

Parameters:
  • tol (float) – tolerance parameter for neighbor identification (default: 0.1).
  • cutoff (float) – cutoff radius in Angstrom to look for trial near-neighbor sites (default: 10.0).
get_nn_info(structure, n)[source]

Get all near-neighbor sites as well as the associated image locations and weights of the site with index n using the closest neighbor distance-based method.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine near neighbors.
Returns:

tuples, each one

of which represents a neighbor site, its image location, and its weight.

Return type:

siw (list of tuples (Site, array, float))

class MinimumOKeeffeNN(tol=0.1, cutoff=10.0)[source]

Bases: pymatgen.analysis.local_env.NearNeighbors

Determine near-neighbor sites and coordination number using the neighbor(s) at closest relative distance, d_min_OKeffee, plus some relative tolerance, where bond valence parameters from O’Keeffe’s bond valence method (J. Am. Chem. Soc. 1991, 3226-3229) are used to calculate relative distances.

Parameters:
  • tol (float) – tolerance parameter for neighbor identification (default: 0.1).
  • cutoff (float) – cutoff radius in Angstrom to look for trial near-neighbor sites (default: 10.0).
get_nn_info(structure, n)[source]

Get all near-neighbor sites as well as the associated image locations and weights of the site with index n using the closest relative neighbor distance-based method with O’Keeffe parameters.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine near neighbors.
Returns:

tuples, each one

of which represents a neighbor site, its image location, and its weight.

Return type:

siw (list of tuples (Site, array, float))

class MinimumVIRENN(tol=0.1, cutoff=10.0)[source]

Bases: pymatgen.analysis.local_env.NearNeighbors

Determine near-neighbor sites and coordination number using the neighbor(s) at closest relative distance, d_min_VIRE, plus some relative tolerance, where atom radii from the ValenceIonicRadiusEvaluator (VIRE) are used to calculate relative distances.

Parameters:
  • tol (float) – tolerance parameter for neighbor identification (default: 0.1).
  • cutoff (float) – cutoff radius in Angstrom to look for trial near-neighbor sites (default: 10.0).
get_nn_info(structure, n)[source]

Get all near-neighbor sites as well as the associated image locations and weights of the site with index n using the closest relative neighbor distance-based method with VIRE atomic/ionic radii.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine near neighbors.
Returns:

tuples, each one

of which represents a neighbor site, its image location, and its weight.

Return type:

siw (list of tuples (Site, array, float))

class NearNeighbors[source]

Bases: object

Base class to determine near neighbors that typically include nearest neighbors and others that are within some tolerable distance.

get_cn(structure, n, use_weights=False)[source]

Get coordination number, CN, of site with index n in structure.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine CN.
  • use_weights (boolean) – flag indicating whether (True) to use weights for computing the coordination number or not (False, default: each coordinated site has equal weight).
Returns:

coordination number.

Return type:

cn (integer or float)

get_nn(structure, n)[source]

Get near neighbors of site with index n in structure.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site in structure for which to determine neighbors.
Returns:

near neighbors.

Return type:

sites (list of Site objects)

get_nn_images(structure, n)[source]

Get image location of all near neighbors of site with index n in structure.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine the image location of near neighbors.
Returns:

image locations of

near neighbors.

Return type:

images (list of 3D integer array)

get_nn_info(structure, n)[source]

Get all near-neighbor sites as well as the associated image locations and weights of the site with index n.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine near-neighbor information.
Returns:

each dictionary provides information

about a single near neighbor, where key ‘site’ gives access to the corresponding Site object, ‘image’ gives the image location, and ‘weight’ provides the weight that a given near-neighbor site contributes to the coordination number (1 or smaller), ‘site_index’ gives index of the corresponding site in the original structure.

Return type:

siw (list of dicts)

get_weights_of_nn_sites(n)[source]

Get weight associated with each near neighbor of site with index n in structure.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine the weights.
Returns:

near-neighbor weights.

Return type:

weights (list of floats)

class ValenceIonicRadiusEvaluator(structure)[source]

Bases: object

Computes site valences and ionic radii for a structure using bond valence analyzer

Parameters:structure – pymatgen.core.structure.Structure
radii

List of ionic radii of elements in the order of sites.

structure

Returns oxidation state decorated structure.

valences

List of oxidation states of elements in the order of sites.

class VoronoiNN(tol=0, targets=None, cutoff=10.0, allow_pathological=False)[source]

Bases: pymatgen.analysis.local_env.NearNeighbors

Uses a Voronoi algorithm to determine near neighbors for each site in a structure.

Parameters:
  • tol (float) – tolerance parameter for near-neighbor finding (default: 0).
  • targets (Element or list of Elements) – target element(s).
  • cutoff (float) – cutoff radius in Angstrom to look for near-neighbor atoms. Defaults to 10.0.
  • allow_pathological (bool) – whether to allow infinite vertices in determination of Voronoi coordination.
get_nn_info(structure, n)[source]

” Get all near-neighbor sites as well as the associated image locations and weights of the site with index n in structure using Voronoi decomposition.

Parameters:
  • structure (Structure) – input structure.
  • n (integer) – index of site for which to determine near-neighbor sites.
Returns:

tuples, each one

of which represents a coordinated site, its image location, and its weight.

Return type:

siw (list of tuples (Site, array, float))

get_voronoi_polyhedra(structure, n)[source]

Gives a weighted polyhedra around a site. This uses the Voronoi construction with solid angle weights. See ref: A Proposed Rigorous Definition of Coordination Number, M. O’Keeffe, Acta Cryst. (1979). A35, 772-775

Parameters:
  • structure (Structure) – structure for which to evaluate the coordination environment.
  • n (integer) – site index.
Returns:

A dict of sites sharing a common Voronoi facet with the site n and their solid angle weights

get_neighbors_of_site_with_index(struct, n, approach='min_dist', delta=0.1, cutoff=10.0)[source]

Returns the neighbors of a given site using a specific neighbor-finding method.

Parameters:
  • struct (Structure) – input structure.
  • n (int) – index of site in Structure object for which motif type is to be determined.
  • approach (str) – type of neighbor-finding approach, where “min_dist” will use the MinimumDistanceNN class, “voronoi” the VoronoiNN class, “min_OKeeffe” the MinimumOKeeffe class, and “min_VIRE” the MinimumVIRENN class.
  • delta (float) – tolerance involved in neighbor finding.
  • cutoff (float) – (large) radius to find tentative neighbors.

Returns: neighbor sites.

get_okeeffe_distance_prediction(el1, el2)[source]

Returns an estimate of the bond valence parameter (bond length) using the derived parameters from ‘Atoms Sizes and Bond Lengths in Molecules and Crystals’ (O’Keeffe & Brese, 1991). The estimate is based on two experimental parameters: r and c. The value for r is based off radius, while c is (usually) the Allred-Rochow electronegativity. Values used are not generated from pymatgen, and are found in ‘okeeffe_params.json’.

Parameters:el2 (el1,) – two Element objects
Returns:a float value of the predicted bond length
get_okeeffe_params(el_symbol)[source]

Returns the elemental parameters related to atom size and electronegativity which are used for estimating bond-valence parameters (bond length) of pairs of atoms on the basis of data provided in ‘Atoms Sizes and Bond Lengths in Molecules and Crystals’ (O’Keeffe & Brese, 1991).

Parameters:el_symbol (str) – element symbol.
Returns:
atom-size (‘r’) and electronegativity-related (‘c’)
parameter.
Return type:(dict)
site_is_of_motif_type(struct, n, approach='min_dist', delta=0.1, cutoff=10.0, thresh=None)[source]

Returns the motif type of the site with index n in structure struct; currently featuring “tetrahedral”, “octahedral”, “bcc”, and “cp” (close-packed: fcc and hcp) as well as “square pyramidal” and “trigonal bipyramidal”. If the site is not recognized, “unrecognized” is returned. If a site should be assigned to two different motifs, “multiple assignments” is returned.

Parameters:
  • struct (Structure) – input structure.
  • n (int) – index of site in Structure object for which motif type is to be determined.
  • approach (str) – type of neighbor-finding approach, where “min_dist” will use the MinimumDistanceNN class, “voronoi” the VoronoiNN class, “min_OKeeffe” the MinimumOKeeffe class, and “min_VIRE” the MinimumVIRENN class.
  • delta (float) – tolerance involved in neighbor finding.
  • cutoff (float) – (large) radius to find tentative neighbors.
  • thresh (dict) – thresholds for motif criteria (currently, required keys and their default values are “qtet”: 0.5, “qoct”: 0.5, “qbcc”: 0.5, “q6”: 0.4).

Returns: motif type (str).

solid_angle(center, coords)[source]

Helper method to calculate the solid angle of a set of coords from the center.

Parameters:
  • center (3x1 array) – Center to measure solid angle from.
  • coords (Nx3 array) – List of coords to determine solid angle.
Returns:

The solid angle.