pymatgen.analysis.diffusion.neb.full_path_mapper module

Migration Graph Analysis.

class ChargeBarrierGraph(structure: Structure, m_graph: StructureGraph, potential_field: VolumetricData, potential_data_key: str, **kwargs)[source]

Bases: MigrationGraph

A Migration graph with additional charge density analysis on the charge density of the host material.

Construct the MigrationGraph object using a VolumetricData object. The graph is constructed using the structure, and cost values are assigned based on charge density analysis.

Parameters:
  • structure (Structure) – Input structure.

  • m_graph (StructureGraph) – Input structure graph.

  • potential_field – Input VolumetricData object that describes the field does not have to contains all the metastable sites.

  • potential_data_key (str) – Key for potential data.

  • **kwargs – Passthru for kwargs.

get_least_chg_path()[source]

obtain an intercolating pathway through the material that has the least amount of charge :returns: list of hops.

get_summary_dict(add_keys: list[str] | None = None)[source]

Dictionary format, for saving to database.

populate_edges_with_chg_density_info(tube_radius=1)[source]
Parameters:

tube_radius – Tube radius.

class MigrationGraph(structure: Structure, m_graph: StructureGraph, symprec=0.1, vac_mode=False)[source]

Bases: MSONable

A python object for handling the migratrion graph of a given base structure and mobile specie sites within that base structure.

Each mobile site is a node in the migration graph and hops (between sites) are edges. The migration graph object uses the Structure Graph from pymatgen.

Construct the MigrationGraph object using a potential_field will all mobile sites occupied. A potential_field graph is generated by connecting all sites as specified by the migration graph. The sites are decorated with Migration graph objects and then grouped together based on their equivalence.

Parameters:
  • structure – Structure with base framework and mobile sites. When used with structure_is_base = True, only the base framework structure, does not contain any migrating sites.

  • m_graph – The StructureGraph object that defines the migration network

  • symprec (float) – Symmetry precision to determine equivalence of migration events

  • vac_mode (Bool) – indicates whether vacancy mode should be used

add_data_to_similar_edges(target_label: int | str, data: dict, m_hop: MigrationHop | None = None)[source]

Insert data to all edges with the same label :param target_label: The edge uniqueness label are adding data :param data: The data to passed to the different edges :param m_hop: If the data is an array, and m_hop is set, it uses the reference migration path to :param determine whether the data needs to be flipped so that 0–>1 is different from 1–>0.:

assign_cost_to_graph(cost_keys=None)[source]

Read the data dict on each add and populate a cost key :param cost_keys: a list of keys for data on each edge.

The SC Graph is decorated with a “cost” key that is the product of the different keys here.

get_path(max_val=100000, flip_hops=True)[source]

obtain a pathway through the material using hops that are in the current graph Basic idea:

Get an endpoint p1 in the graph that is outside the current unit cell Ask the graph for a pathway that connects to p1 from either within the (0,0,0) cell or any other neighboring UC not containing p1.

Parameters:
  • max_val – Filter the graph by a cost

  • flip_hops – If true, hops in paths returned will be flipped so isites and esites match to form a coherent path. If false, hops will retain their original orientation from the migration graph.

Returns:

Each dict contains the information of a hop

Return type:

Generator for list of Dicts

static get_structure_from_entries(entries: list[ComputedStructureEntry], migrating_ion_entry: ComputedEntry, **kwargs) list[Structure][source]

Read in a list of base entries and inserted entries. Return a list of structures that contains metastable sites for the migration species decorated with a “insertion_energy” property.

Parameters:
  • entries – list of entries, must contain a mixture of inserted and empty structures.

  • migrating_ion_entry – The metallic phase of the working ion, used to calculate insertion energies.

  • **kwargs – Passthrough for kwargs.

Additional Kwargs:

symprec: symmetry parameter for SpacegroupAnalyzer ltol: Fractional length tolerance for StructureMatcher stol: Site tolerance for StructureMatcher angle_tol: Angle tolerance for StructureMatcher and SpacegroupAnalyzer only_single_cat: If True, only use single cation insertions so the site energy is more accurate use_strict_tol: halve the ltol and stol parameter for more strict matching.

Returns:

a list of host structures with working ion on all the metastable sites. The structures are ranked by the number of metastable sites (most is first) If the given entries are not enough to construct such a structure, return an empty list.

get_summary_dict(added_keys: list[str] | None = None) dict[source]

Dictionary format, for saving to database.

Parameters:

added_keys – a list of keys for data on each edge.

Returns:

Dict.

property host_structure: Structure

A structure that only contains the non-migrating species.

property only_sites: Structure

A structure that only contains the migrating species.

property symm_structure: SymmetrizedStructure

The symmetrized structure with the present item’s symprec value.

property unique_hops

The unique hops dictionary keyed by the hop label.

classmethod with_base_structure(base_structure: Structure, m_graph: StructureGraph, **kwargs) MigrationGraph[source]
Parameters:
  • base_structure – base framework structure that does not contain any migrating sites.

  • m_graph – The StructureGraph object that defines the migration network.

  • **kwargs – Passthrough for kwargs.

Returns:

A constructed MigrationGraph object

classmethod with_distance(structure: Structure, migrating_specie: str, max_distance: float, **kwargs) MigrationGraph[source]

Using a specific nn strategy to get the connectivity graph between all the migrating ion sites.

Parameters:
  • structure – Input structure that contains all sites.

  • migrating_specie – The specie that migrates. E.g. “Li”.

  • max_distance – Maximum length of NEB path in the unit of Angstrom. Defaults to None, which means you are setting the value to the min cutoff until finding 1D or >1D percolating paths.

  • **kwargs – Passthrough for kwargs.

Returns:

A constructed MigrationGraph object

classmethod with_local_env_strategy(structure: Structure, migrating_specie: str, nn: NearNeighbors, **kwargs) MigrationGraph[source]

Using a specific nn strategy to get the connectivity graph between all the migrating ion sites.

Parameters:
  • structure – Input structure that contains all sites.

  • migrating_specie – The specie that migrates. E.g. “Li”.

  • nn – The specific local environment object used to connect the migrating ion sites.

  • **kwargs – Passthrough for kwargs.

Returns:

A constructed MigrationGraph object

almost(a, b)[source]

Return true if the values are almost equal.

check_uc_hop(sc_hop, uc_hop)[source]

See if hop in the 2X2X2 supercell and a unit cell hop are equivalent under lattice translation.

Parameters:
  • sc_hop – MigrationHop object form pymatgen-diffusion.

  • uc_hop – MigrationHop object form pymatgen-diffusion.

Returns:

image vector of length 3 Is the UC hop flip of the SC hop

generic_groupby(list_in: list, comp: ~typing.Callable = <built-in function eq>)[source]

Group a list of unsortable objects.

Parameters:
  • list_in – A list of generic objects

  • comp – (Default value = operator.eq) The comparator

Returns:

[int] list of labels for the input list

get_hop_site_sequence(hop_list: list[dict], start_u: int | str, key: str | None = None) list[source]

Read in a list of hop dictionaries and print the sequence of sites (and relevant property values if any).

Parameters:
  • hop_list – a list of the data on a sequence of hops

  • start_u – the site index of the starting sites

  • key (optional) – property to track in a hop (e.g.: “hop_distance”)

Returns:

String representation of the hop sequence (and property values if any)

get_only_sites_from_structure(structure: Structure, migrating_specie: str) Structure[source]

Get a copy of the structure with only the migrating sites.

Parameters:
  • structure – The full_structure that contains all the sites

  • migrating_specie – The name of migrating species

Returns:

Structure with all possible migrating ion sites

Return type:

Structure

map_hop_sc2uc(sc_hop: MigrationHop, mg: MigrationGraph)[source]

Map a given hop in the SC onto the UC.

Parameters:
  • sc_hop – MigrationHop object form pymatgen-diffusion.

  • mg – MigrationGraph object from pymatgen-diffusion.

Note

For now assume that the SC is exactly 2x2x2 of the UC. Can add in the parsing of different SC’s later For a migration event in the SC from (0.45,0,0)–>(0.55,0,0) the UC hop might be (0.9,0,0)–>(0.1,0,0)[img:100] for the inverse of (0.1,0,0)–>(-0.1,0,0) the code needs to account for both those cases

order_path(hop_list: list[dict], start_u: int | str) list[dict][source]

Takes a list of hop dictionaries and flips hops (switches isite and esite) as needed to form a coherent path / sequence of sites according to get_hop_site_sequence(). For example if hop_list = [{iindex:0, eindex:1, etc.}, {iindex:0, eindex:1, etc.}] then the output is [{iindex:0, eindex:1, etc.}, {iindex:1, eindex:0, etc.}] so that the following hop iindex matches the previous hop’s eindex.

Parameters:
  • hop_list – a list of the data on a sequence of hops

  • start_u – the site index of the starting sites

Returns:

a list of the data on a sequence of hops with hops in coherent orientation