pymatgen.analysis.structure_analyzer module

class JMolCoordFinder(el_radius_updates=None)[source]

Bases: object

Determine coordinated 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.

Initialize coordination finder parameters (atomic radii)

Parameters:el_radius_updates – (dict) symbol->float to override default atomic radii table values
get_coordinated_sites(structure, n, tol=0.001)[source]

Get the coordinated sites for a site :param structure: (Structure) :param n: (int) index of site in the structure to analyze :param tol: (float) a numerical tolerance to extend search

Returns: ([sites]) a list of coordinated sites

get_coordination_number(structure, n, tol=0.001)[source]

Get the coordination number of a site :param structure: (Structure) :param n: (int) index of site in the structure to get CN for :param tol: (float) a numerical tolerance to extend search

Returns: (int) the coordination number

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

class OrderParameters(types, parameters=None, cutoff=-10.0)[source]

Bases: object

This class permits the calculation of various types of local structure order parameters.

Parameters:
  • types ([string]) –

    list of strings representing the types of order parameters to be calculated. Note that multiple mentions of the same type may occur. Currently available types recognize following environments:

    ”cn”: simple coordination number—normalized
    if desired;

    ”sgl_bd”: single bonds; “bent”: bent (angular) coordinations

    (Zimmermann & Jain, in progress, 2017);

    ”T”: T-shape coordinations; “see_saw”: see saw-like coordinations; “tet”: tetrahedra

    (Zimmermann et al., submitted, 2017);
    ”oct”: octahedra
    (Zimmermann et al., submitted, 2017);
    ”bcc”: body-centered cubic environments (Peters,
    1. Chem. Phys., 131, 244103, 2009);

    ”tri_plan”: trigonal planar environments; “sq_plan”: square planar environments; “pent_plan”: pentagonal planar environments; “tri_pyr”: trigonal pyramids (coordinated atom is in

    the center of the basal plane);

    ”sq_pyr”: square pyramids; “pent_pyr”: pentagonal pyramids; “hex_pyr”: hexagonal pyramids; “tri_bipyr”: trigonal bipyramids; “sq_bipyr”: square bipyramids; “pent_bipyr”: pentagonal bipyramids; “hex_bipyr”: hexagonal bipyramids; “cuboct”: cuboctahedra; “q2”: motif-unspecific bond orientational order

    parameter (BOOP) of weight l=2 (Steinhardt et al., Phys. Rev. B, 28, 784-805, 1983);

    ”q4”: BOOP of weight l=4; “q6”: BOOP of weight l=6. “reg_tri”: regular triangle with varying height

    to basal plane;

    ”sq”: square coordination (cf., “reg_tri”); “oct_legacy”: original Peters-style OP recognizing

    octahedral coordination environments (Zimmermann et al., J. Am. Chem. Soc., 137, 13352-13361, 2015) that can, however, produce small negative values sometimes.

    ”sq_pyr_legacy”: square pyramids (legacy);

  • parameters ([dict]) –

    list of dictionaries that store float-type parameters associated with the definitions of the different order parameters (length of list = number of OPs). If an entry is None, default values are used that are read from the op_paras.yaml file. With few exceptions, 9 different parameters are used across all OPs:

    ”norm”: normalizing constant (used in “cn”
    (default value: 1)).
    ”TA”: target angle (TA) in fraction of 180 degrees
    (“bent” (1), “tet” (0.6081734479693927), “tri_plan” (0.66666666667), “pent_plan” (0.6), “sq_pyr_legacy” (0.5)).
    ”IGW_TA”: inverse Gaussian width (IGW) for penalizing
    angles away from the target angle in inverse fractions of 180 degrees to (“bent” and “tet” (15), “tri_plan” (13.5), “pent_plan” (18), “sq_pyr_legacy” (30)).
    ”IGW_EP”: IGW for penalizing angles away from the
    equatorial plane (EP) at 90 degrees (“T”, “see_saw”, “oct”, “sq_plan”, “tri_pyr”, “sq_pyr”, “pent_pyr”, “hex_pyr”, “tri_bipyr”, “sq_bipyr”, “pent_bipyr”, “hex_bipyr”, and “oct_legacy” (18)).
    ”fac_AA”: factor applied to azimuth angle (AA) in cosine
    term (“T”, “tri_plan”, and “sq_plan” (1), “tet”, “tri_pyr”, and “tri_bipyr” (1.5), “oct”, “sq_pyr”, “sq_bipyr”, and “oct_legacy” (2), “pent_pyr” and “pent_bipyr” (2.5), “hex_pyr” and “hex_bipyr” (3)).
    ”exp_cos_AA”: exponent applied to cosine term of AA
    (“T”, “tet”, “oct”, “tri_plan”, “sq_plan”, “tri_pyr”, “sq_pyr”, “pent_pyr”, “hex_pyr”, “tri_bipyr”, “sq_bipyr”, “pent_bipyr”, “hex_bipyr”, and “oct_legacy” (2)).
    ”min_SPP”: smallest angle (in radians) to consider
    a neighbor to be at South pole position (“see_saw”, “oct”, “bcc”, “sq_plan”, “tri_bipyr”, “sq_bipyr”, “pent_bipyr”, “hex_bipyr”, “cuboct”, and “oct_legacy” (2.792526803190927)).
    ”IGW_SPP”: IGW for penalizing angles away from South
    pole position (“see_saw”, “oct”, “bcc”, “sq_plan”, “tri_bipyr”, “sq_bipyr”, “pent_bipyr”, “hex_bipyr”, “cuboct”, and “oct_legacy” (15)).
    ”w_SPP”: weight for South pole position relative to
    equatorial positions (“see_saw” and “sq_plan” (1), “cuboct” (1.8), “tri_bipyr” (2), “oct”, “sq_bipyr”, and “oct_legacy” (3), “pent_bipyr” (4), “hex_bipyr” (5), “bcc” (6)).
  • cutoff (float) – Cutoff radius to determine which nearest neighbors are supposed to contribute to the order parameters. If the value is negative the neighboring sites found by distance and cutoff radius are further pruned using the get_coordinated_sites method from the VoronoiCoordFinder class.
compute_trigonometric_terms(thetas, phis)[source]

” Computes trigonometric terms that are required to calculate bond orientational order parameters using internal variables.

Parameters:
  • thetas ([float]) – polar angles of all neighbors in radians.
  • phis ([float]) – azimuth angles of all neighbors in radians. The list of azimuth angles of all neighbors in radians. The list of azimuth angles is expected to have the same size as the list of polar angles; otherwise, a ValueError is raised. Also, the two lists of angles have to be coherent in order. That is, it is expected that the order in the list of azimuth angles corresponds to a distinct sequence of neighbors. And, this sequence has to equal the sequence of neighbors in the list of polar angles.
get_order_parameters(structure, n, indices_neighs=None, tol=0.0, target_spec=None)[source]

Compute all order parameters of site n.

Parameters:
  • structure (Structure) – input structure.
  • n (int) – index of site in input structure, for which OPs are to be calculated. Note that we do not use the sites iterator here, but directly access sites via struct[index].
  • indices_neighs ([int]) – list of indices of those neighbors in Structure object structure that are to be considered for OP computation. This optional argument overwrites the way neighbors are to be determined as defined in the constructor (i.e., Voronoi coordination finder via negative cutoff radius vs constant cutoff radius if cutoff was positive). We do not use information about the underlying structure lattice if the neighbor indices are explicitly provided. This has two important consequences. First, the input Structure object can, in fact, be a simple list of Site objects. Second, no nearest images of neighbors are determined when providing an index list. Note furthermore that this neighbor determination type ignores the optional target_spec argument.
  • tol (float) – threshold of weight (= solid angle / maximal solid angle) to determine if a particular pair is considered neighbors; this is relevant only in the case when Voronoi polyhedra are used to determine coordination
  • target_spec (Specie) – target species to be considered when calculating the order parameters of site n; None includes all species of input structure.
Returns:

representing order parameters. Should it not be possible to compute a given OP for a conceptual reason, the corresponding entry is None instead of a float. For Steinhardt et al.’s bond orientational OPs and the other geometric OPs (“tet”, “oct”, “bcc”, etc.), this can happen if there is a single neighbor around site n in the structure because that does not permit calculation of angles between multiple neighbors.

Return type:

[floats]

get_parameters(index)[source]

Returns list of floats that represents the parameters associated with calculation of the order parameter that was defined at the index provided. Attention: the parameters do not need to equal those originally inputted because of processing out of efficiency reasons.

Parameters:index (int) – index of order parameter for which associated parameters are to be returned.
Returns:parameters of a given OP.
Return type:[float]
get_q2(thetas=None, phis=None)[source]

Calculates the value of the bond orientational order parameter of weight l=2. If the function is called with non-empty lists of polar and azimuthal angles the corresponding trigonometric terms are computed afresh. Otherwise, it is expected that the compute_trigonometric_terms function has been just called.

Parameters:
  • thetas ([float]) – polar angles of all neighbors in radians.
  • phis ([float]) – azimuth angles of all neighbors in radians.
Returns:

bond orientational order parameter of weight l=2

corresponding to the input angles thetas and phis.

Return type:

float

get_q4(thetas=None, phis=None)[source]

Calculates the value of the bond orientational order parameter of weight l=4. If the function is called with non-empty lists of polar and azimuthal angles the corresponding trigonometric terms are computed afresh. Otherwise, it is expected that the compute_trigonometric_terms function has been just called.

Parameters:
  • thetas ([float]) – polar angles of all neighbors in radians.
  • phis ([float]) – azimuth angles of all neighbors in radians.
Returns:

bond orientational order parameter of weight l=4

corresponding to the input angles thetas and phis.

Return type:

float

get_q6(thetas=None, phis=None)[source]

Calculates the value of the bond orientational order parameter of weight l=6. If the function is called with non-empty lists of polar and azimuthal angles the corresponding trigonometric terms are computed afresh. Otherwise, it is expected that the compute_trigonometric_terms function has been just called.

Parameters:
  • thetas ([float]) – polar angles of all neighbors in radians.
  • phis ([float]) – azimuth angles of all neighbors in radians.
Returns:

bond orientational order parameter of weight l=6

corresponding to the input angles thetas and phis.

Return type:

float

get_type(index)[source]

Return type of order parameter at the index provided and represented by a short string.

Parameters:index (int) – index of order parameter for which type is to be returned.
Returns:OP type.
Return type:str
last_nneigh

” :returns:

the number of neighbors encountered during the most
recent order parameter calculation. A value of -1 indicates that no such calculation has yet been performed for this instance.
Return type:int
num_ops

” :returns:

the number of different order parameters that are targeted
to be calculated.
Return type:int
class OxideType(structure, relative_cutoff=1.1)[source]

Bases: object

Separate class for determining oxide type.

Parameters:
  • structure – Input structure.
  • relative_cutoff – Relative_cutoff * act. cutoff stipulates the max. distance two O atoms must be from each other. Default value is 1.1. At most 1.1 is recommended, nothing larger, otherwise the script cannot distinguish between superoxides and peroxides.
parse_oxide()[source]

Determines if an oxide is a peroxide/superoxide/ozonide/normal oxide.

Returns:Type of oxide ozonide/peroxide/superoxide/hydroxide/None. nbonds (int): Number of peroxide/superoxide/hydroxide bonds in structure.
Return type:oxide_type (str)
class RelaxationAnalyzer(initial_structure, final_structure)[source]

Bases: object

This class analyzes the relaxation in a calculation.

Please note that the input and final structures should have the same ordering of sites. This is typically the case for most computational codes.

Parameters:
  • initial_structure (Structure) – Initial input structure to calculation.
  • final_structure (Structure) – Final output structure from calculation.
get_percentage_bond_dist_changes(max_radius=3.0)[source]

Returns the percentage bond distance changes for each site up to a maximum radius for nearest neighbors.

Parameters:max_radius (float) – Maximum radius to search for nearest neighbors. This radius is applied to the initial structure, not the final structure.
Returns:Bond distance changes as a dict of dicts. E.g., {index1: {index2: 0.011, …}}. For economy of representation, the index1 is always less than index2, i.e., since bonding between site1 and siten is the same as bonding between siten and site1, there is no reason to duplicate the information or computation.
get_percentage_lattice_parameter_changes()[source]

Returns the percentage lattice parameter changes.

Returns:A dict of the percentage change in lattice parameter, e.g., {‘a’: 0.012, ‘b’: 0.021, ‘c’: -0.031} implies a change of 1.2%, 2.1% and -3.1% in the a, b and c lattice parameters respectively.
get_percentage_volume_change()[source]

Returns the percentage volume change.

Returns:Volume change in percentage, e.g., 0.055 implies a 5.5% increase.
class VoronoiAnalyzer(cutoff=5.0, qhull_options='Qbb Qc Qz')[source]

Bases: object

Performs a statistical analysis of Voronoi polyhedra around each site. Each Voronoi polyhedron is described using Schaefli notation. That is a set of indices {c_i} where c_i is the number of faces with i number of vertices. E.g. for a bcc crystal, there is only one polyhedron notation of which is [0,6,0,8,0,0,…]. In perfect crystals, these also corresponds to the Wigner-Seitz cells. For distorted-crystals, liquids or amorphous structures, rather than one-type, there is a statistical distribution of polyhedra. See ref: Microstructure and its relaxation in Fe-B amorphous system simulated by molecular dynamics,

Stepanyuk et al., J. Non-cryst. Solids (1993), 159, 80-87.
Parameters:
  • cutoff (float) – cutoff distance to search for neighbors of a given atom (default = 5.0)
  • qhull_options (str) – options to pass to qhull (optional)
analyze(structure, n=0)[source]

Performs Voronoi analysis and returns the polyhedra around atom n in Schlaefli notation.

Parameters:
  • structure (Structure) – structure to analyze
  • n (int) – index of the center atom in structure
Returns:

<c3,c4,c6,c6,c7,c8,c9,c10>

where c_i denotes number of facets with i vertices.

Return type:

voronoi index of n

analyze_structures(structures, step_freq=10, most_frequent_polyhedra=15)[source]

Perform Voronoi analysis on a list of Structures. Note that this might take a significant amount of time depending on the size and number of structures.

Parameters:
  • structures (list) – list of Structures
  • (float (cutoff) – cutoff distance around an atom to search for neighbors
  • step_freq (int) – perform analysis every step_freq steps
  • qhull_options (str) – options to pass to qhull
  • most_frequent_polyhedra (int) – this many unique polyhedra with highest frequences is stored.
Returns:

A list of tuples in the form (voronoi_index,frequency)

static plot_vor_analysis(voronoi_ensemble)[source]
class VoronoiConnectivity(structure, cutoff=10)[source]

Bases: object

Computes the solid angles swept out by the shared face of the voronoi polyhedron between two sites.

Parameters:
  • structure (Structure) – Input structure
  • cutoff (float) –
connectivity_array

Provides connectivity array.

Returns:An array of shape [atomi, atomj, imagej]. atomi is the index of the atom in the input structure. Since the second atom can be outside of the unit cell, it must be described by both an atom index and an image index. Array data is the solid angle of polygon between atomi and imagej of atomj
Return type:connectivity
get_connections()[source]

Returns a list of site pairs that are Voronoi Neighbors, along with their real-space distances.

get_sitej(site_index, image_index)[source]

Assuming there is some value in the connectivity array at indices (1, 3, 12). sitei can be obtained directly from the input structure (structure[1]). sitej can be obtained by passing 3, 12 to this function

Parameters:
  • site_index (int) – index of the site (3 in the example)
  • image_index (int) – index of the image (12 in the example)
max_connectivity

returns the 2d array [sitei, sitej] that represents the maximum connectivity of site i to any periodic image of site j

class VoronoiCoordFinder(structure, target=None, cutoff=10.0, allow_pathological=False)[source]

Bases: object

Uses a Voronoi algorithm to determine the coordination for each site in a structure.

Parameters:
  • structure (Structure) – Input structure
  • target ([Element/Specie]) – A list of target species to determine coordination for.
  • cutoff (float) – Radius in Angstrom cutoff to look for coordinating atoms. Defaults to 10.0.
  • allow_pathological (bool) – whether to allow infinite vertices in determination of Voronoi coordination
get_coordinated_sites(n, tol=0, target=None)[source]

Returns the sites that are in the coordination radius of site with index n.

Parameters:
  • n (int) – Site index.
  • tol (float) – Weight tolerance to determine if a particular pair is considered a neighbor.
  • target (Element) – Target element
Returns:

Sites coordinating input site.

get_coordination_number(n)[source]

Returns the coordination number of site with index n.

Parameters:n (int) – Site index
get_voronoi_polyhedra(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:n (int) – Site index
Returns:A dict of sites sharing a common Voronoi facet with the site n and their solid angle weights
average_coordination_number(structures, freq=10)[source]

Calculates the ensemble averaged Voronoi coordination numbers of a list of Structures using VoronoiCoordFinder. Typically used for analyzing the output of a Molecular Dynamics run. :param structures: list of Structures. :type structures: list :param freq: sampling frequency of coordination number [every freq steps]. :type freq: int

Returns:Dictionary of elements as keys and average coordination numbers as values.
contains_peroxide(structure, relative_cutoff=1.1)[source]

Determines if a structure contains peroxide anions.

Parameters:
  • structure (Structure) – Input structure.
  • relative_cutoff – The peroxide bond distance is 1.49 Angstrom. Relative_cutoff * 1.49 stipulates the maximum distance two O atoms must be to each other to be considered a peroxide.
Returns:

Boolean indicating if structure contains a peroxide anion.

get_dimensionality(structure, max_hkl=2, el_radius_updates=None, min_slab_size=5, min_vacuum_size=5, standardize=True, bonds=None)[source]

This method returns whether a structure is 3D, 2D (layered), or 1D (linear chains or molecules) according to the algorithm published in Gorai, P., Toberer, E. & Stevanovic, V. Computational Identification of Promising Thermoelectric Materials Among Known Quasi-2D Binary Compounds. J. Mater. Chem. A 2, 4136 (2016).

Note that a 1D structure detection might indicate problems in the bonding algorithm, particularly for ionic crystals (e.g., NaCl)

Users can change the behavior of bonds detection by passing either el_radius_updates to update atomic radii for auto-detection of max bond distances, or bonds to explicitly specify max bond distances for atom pairs. Note that if you pass both, el_radius_updates are ignored.

Parameters:
  • structure – (Structure) structure to analyze dimensionality for
  • max_hkl – (int) max index of planes to look for layers
  • el_radius_updates – (dict) symbol->float to update atomic radii
  • min_slab_size – (float) internal surface construction parameter
  • min_vacuum_size – (float) internal surface construction parameter
  • standardize (bool) – whether to standardize the structure before analysis. Set to False only if you already have the structure in a convention where layers / chains will be along low <hkl> indexes.
  • bonds ({(specie1, specie2) – max_bond_dist}: bonds are specified as a dict of tuples: float of specie1, specie2 and the max bonding distance. For example, PO4 groups may be defined as {(“P”, “O”): 3}.
Returns: (int) the dimensionality of the structure - 1 (molecules/chains),
2 (layered), or 3 (3D)
get_max_bond_lengths(structure, el_radius_updates=None)[source]

Provides max bond length estimates for a structure based on the JMol table and algorithms.

Parameters:
  • structure – (structure)
  • el_radius_updates – (dict) symbol->float to update atomic radii
Returns: (dict) - (Element1, Element2) -> float. The two elements are
ordered by Z.
gramschmidt(vin, uin)[source]

Returns that part of the first input vector that is orthogonal to the second input vector. The output vector is not normalized.

Parameters:
  • vin (numpy array) – first input vector
  • uin (numpy array) – second input vector
oxide_type(structure, relative_cutoff=1.1, return_nbonds=False)[source]

Determines if an oxide is a peroxide/superoxide/ozonide/normal oxide

Parameters:
  • structure (Structure) – Input structure.
  • relative_cutoff (float) – Relative_cutoff * act. cutoff stipulates the max distance two O atoms must be from each other.
  • return_nbonds (bool) – Should number of bonds be requested?
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.

sulfide_type(structure)[source]

Determines if a structure is a sulfide/polysulfide

Parameters:structure (Structure) – Input structure.
Returns:(str) sulfide/polysulfide/sulfate