pymatgen.analysis.magnetism package
Package for analysis of magnetic structures.
Submodules
pymatgen.analysis.magnetism.analyzer module
This module provides some useful functions for dealing with magnetic Structures (e.g. Structures with associated magmom tags).
- class CollinearMagneticStructureAnalyzer(structure: Structure, overwrite_magmom_mode: str | OverwriteMagmomMode = OverwriteMagmomMode.none, round_magmoms: bool = False, detect_valences: bool = False, make_primitive: bool = True, default_magmoms: dict | None = None, set_net_positive: bool = True, threshold: float = 0, threshold_nonmag: float = 0.1, threshold_ordering: float = 1e-08)[source]
Bases:
object
A class which provides a few helpful methods to analyze collinear magnetic structures.
If magnetic moments are not defined, moments will be taken either from default_magmoms.yaml (similar to the default magmoms in MPRelaxSet, with a few extra definitions) or from a specie:magmom dict provided by the default_magmoms kwarg.
Input magmoms can be replaced using the ‘overwrite_magmom_mode’ kwarg. This can be: * “none” to do nothing, * “respect_sign” which will overwrite existing magmoms with
those from default_magmoms but will keep sites with positive magmoms positive, negative magmoms negative and zero magmoms zero,
“respect_zeros”, which will give a ferromagnetic structure (all positive magmoms from default_magmoms) but still keep sites with zero magmoms as zero,
“replace_all” which will try to guess initial magmoms for all sites in the structure irrespective of input structure (this is most suitable for an initial DFT calculation),
“replace_all_if_undefined” is the same as “replace_all” but only if no magmoms are defined in input structure, otherwise it will respect existing magmoms.
“normalize” will normalize magmoms to unity, but will respect sign (used for comparing orderings), magmoms < threshold will be set to zero
- Parameters:
structure – input Structure object
overwrite_magmom_mode – “respect_sign”, “respect_zeros”, “replace_all”, “replace_all_if_undefined”, “normalize” (default “none”)
round_magmoms – will round input magmoms to specified number of decimal places if integer is supplied, if set to a float will try and group magmoms together using a kernel density estimator of provided width, and extracting peaks of the estimator detect_valences: if True, will attempt to assign valences to input structure
make_primitive – if True, will transform to primitive magnetic cell
default_magmoms – (optional) dict specifying default magmoms
set_net_positive – if True, will change sign of magnetic moments such that the net magnetization is positive. Argument will be ignored if mode “respect_sign” is used.
threshold – number (in Bohr magnetons) below which magmoms will be rounded to zero
threshold_nonmag – number (in Bohr magneton) below which nonmagnetic ions (with no magmom specified in default_magmoms) will be rounded to zero
threshold_ordering – number (absolute of sum of all magmoms, in Bohr magneton) below which total magnetization is treated as zero when defining magnetic ordering. Defaults to 1e-8.
- get_exchange_group_info(symprec: float = 0.01, angle_tolerance: float = 5) tuple[str, int] [source]
Get the information on the symmetry of the Hamiltonian describing the exchange energy of the system, taking into account relative direction of magnetic moments but not their absolute direction.
This is not strictly accurate (e.g. some/many atoms will have zero magnetic moments), but defining symmetry this way is a useful way of keeping track of distinct magnetic orderings within pymatgen.
- Parameters:
symprec – same as SpacegroupAnalyzer (Default value = 1e-2)
angle_tolerance – same as SpacegroupAnalyzer (Default value = 5)
- Returns:
spacegroup_symbol, international_number
- get_ferromagnetic_structure(make_primitive: bool = True) Structure [source]
Get a Structure with all magnetic moments positive or zero.
- Parameters:
make_primitive – Whether to make structure primitive after making all magnetic moments positive (Default value = True)
- Returns:
Structure
- get_nonmagnetic_structure(make_primitive: bool = True) Structure [source]
Get a Structure without magnetic moments defined.
- Parameters:
make_primitive – Whether to make structure primitive after removing magnetic information (Default value = True)
- Returns:
Structure
- get_structure_with_only_magnetic_atoms(make_primitive: bool = True) Structure [source]
Get a Structure with only magnetic atoms present.
- Parameters:
make_primitive – Whether to make structure primitive after removing non-magnetic atoms (Default value = True)
- Returns:
Structure
- get_structure_with_spin() Structure [source]
Get a Structure with species decorated with spin values instead of using magmom site properties.
- property is_magnetic: bool[source]
Convenience property, returns True if any non-zero magmoms present.
- property magnetic_species_and_magmoms: dict[str, Any][source]
A dict of magnetic species and the magnitude of their associated magmoms. Will return a list if there are multiple magmoms per species.
- Returns:
dict of magnetic species and magmoms
- matches_ordering(other: Structure) bool [source]
Compares the magnetic orderings of one structure with another.
- Parameters:
other – Structure to compare
- Returns:
True if magnetic orderings match, False otherwise
- Return type:
bool
- number_of_unique_magnetic_sites(symprec: float = 0.001, angle_tolerance: float = 5) int [source]
- Parameters:
symprec – same as in SpacegroupAnalyzer (Default value = 1e-3)
angle_tolerance – same as in SpacegroupAnalyzer (Default value = 5).
- Returns:
Number of symmetrically-distinct magnetic sites present in structure.
- Return type:
int
- property ordering: Ordering[source]
Apply heuristics to return a magnetic ordering for a collinear magnetic structure. Result is not guaranteed to be correct, just a best guess. Tolerance for minimum total magnetization to be considered ferro/ferrimagnetic is self.threshold_ordering and defaults to 1e-8.
- Returns:
- Enum with values FM: ferromagnetic, FiM: ferrimagnetic,
AFM: antiferromagnetic, NM: non-magnetic or Unknown. Unknown is returned if magnetic moments are not defined or structure is not collinear (in which case a warning is issued).
- Return type:
- property types_of_magnetic_specie: tuple[Element | Species | DummySpecies, ...][source]
Specie->Species rename. Used to maintain backwards compatibility.
- property types_of_magnetic_species: tuple[Element | Species | DummySpecies, ...][source]
Equivalent to Structure.types_of_specie but only returns magnetic species.
- Returns:
types of Species
- Return type:
tuple
- class MagneticDeformation(deformation, type)[source]
Bases:
NamedTuple
Create new instance of MagneticDeformation(deformation, type)
- class MagneticStructureEnumerator(structure: Structure, default_magmoms: dict[str, float] | None = None, strategies: list[str] | tuple[str, ...] = ('ferromagnetic', 'antiferromagnetic'), automatic: bool = True, truncate_by_symmetry: bool = True, transformation_kwargs: dict | None = None)[source]
Bases:
object
Combines MagneticStructureAnalyzer and MagOrderingTransformation to automatically generate a set of transformations for a given structure and produce a list of plausible magnetic orderings.
Generate different collinear magnetic orderings for a given input structure.
If the input structure has magnetic moments defined, it is possible to use these as a hint as to which elements are magnetic, otherwise magnetic elements will be guessed (this can be changed using default_magmoms kwarg).
- Parameters:
structure – input structure
default_magmoms – (optional, defaults provided) dict of magnetic elements to their initial magnetic moments in μB, generally these are chosen to be high-spin since they can relax to a low-spin configuration during a DFT electronic configuration
strategies – different ordering strategies to use, choose from: ferromagnetic, antiferromagnetic, antiferromagnetic_by_motif, ferrimagnetic_by_motif and ferrimagnetic_by_species (here, “motif”, means to use a different ordering parameter for symmetry inequivalent sites)
automatic – if True, will automatically choose sensible strategies
truncate_by_symmetry – if True, will remove very unsymmetrical orderings that are likely physically implausible
transformation_kwargs – keyword arguments to pass to MagOrderingTransformation, to change automatic cell size limits, etc.
- class Ordering(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
Enumeration defining possible magnetic orderings.
- class OverwriteMagmomMode(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
Enumeration defining different modes for analyzer.
- magnetic_deformation(structure_A: Structure, structure_B: Structure) MagneticDeformation [source]
Calculate ‘magnetic deformation proxy’, a measure of deformation (norm of finite strain) between ‘non-magnetic’ (non-spin-polarized) and ferromagnetic structures.
Adapted from Bocarsly et al. 2017, doi: 10.1021/acs.chemmater.6b04729
- Parameters:
structure_A – Structure
structure_B – Structure
- Returns:
MagneticDeformation
pymatgen.analysis.magnetism.heisenberg module
This module implements a simple algorithm for extracting nearest neighbor exchange parameters by mapping low energy magnetic orderings to a Heisenberg model.
- class HeisenbergMapper(ordered_structures, energies, cutoff=0, tol: float = 0.02)[source]
Bases:
object
Compute exchange parameters from low energy magnetic orderings.
Exchange parameters are computed by mapping to a classical Heisenberg model. Strategy is the scheme for generating neighbors. Currently only MinimumDistanceNN is implemented. n+1 unique orderings are required to compute n exchange parameters.
First run a MagneticOrderingsWF to obtain low energy collinear magnetic orderings and find the magnetic ground state. Then enumerate magnetic states with the ground state as the input structure, find the subset of supercells that map to the ground state, and do static calculations for these orderings.
- Parameters:
ordered_structures (list) – Structure objects with magmoms.
energies (list) – Total energies of each relaxed magnetic structure.
cutoff (float) – Cutoff in Angstrom for nearest neighbor search. Defaults to 0 (only NN, no NNN, etc.)
tol (float) – Tolerance (in Angstrom) on nearest neighbor distances being equal.
- estimate_exchange(fm_struct=None, afm_struct=None, fm_e=None, afm_e=None)[source]
Estimate <J> for a structure based on low energy FM and AFM orderings.
- get_exchange()[source]
Take Heisenberg Hamiltonian and corresponding energy for each row and solve for the exchange parameters.
- Returns:
Exchange parameters (meV/atom).
- Return type:
dict[str, float]
- get_heisenberg_model()[source]
Save results of mapping to a HeisenbergModel object.
- Returns:
MSONable object.
- Return type:
- get_interaction_graph(filename=None)[source]
Get a StructureGraph with edges and weights that correspond to exchange interactions and J_ij values, respectively.
- Parameters:
filename (str) – if not None, save interaction graph to filename.
- Returns:
Exchange interaction graph.
- Return type:
- get_low_energy_orderings()[source]
Find lowest energy FM and AFM orderings to compute E_AFM - E_FM.
- Returns:
fm structure with ‘magmom’ site property afm_struct (Structure): afm structure with ‘magmom’ site property fm_e (float): fm energy afm_e (float): afm energy
- Return type:
fm_struct (Structure)
- get_mft_temperature(j_avg)[source]
Crude mean field estimate of critical temperature based on <J> for one sublattice, or solving the coupled equations for a multi-sublattice material.
- Parameters:
j_avg (float) – j_avg (float): Average exchange parameter (meV/atom)
- Returns:
Critical temperature mft_t (K)
- Return type:
float
- class HeisenbergModel(formula=None, structures=None, energies=None, cutoff=None, tol=None, sgraphs=None, unique_site_ids=None, wyckoff_ids=None, nn_interactions=None, dists=None, ex_mat=None, ex_params=None, javg=None, igraph=None)[source]
Bases:
MSONable
Store a Heisenberg model fit to low-energy magnetic orderings. Intended to be generated by HeisenbergMapper.get_heisenberg_model().
- Parameters:
formula (str) – Reduced formula of compound.
structures (list) – Structure objects with magmoms.
energies (list) – Energies of each relaxed magnetic structure.
cutoff (float) – Cutoff in Angstrom for nearest neighbor search.
tol (float) – Tolerance (in Angstrom) on nearest neighbor distances being equal.
sgraphs (list) – StructureGraph objects.
unique_site_ids (dict) – Maps each site to its unique numerical identifier.
wyckoff_ids (dict) – Maps unique numerical identifier to wyckoff position.
nn_interactions (dict) – {i: j} pairs of NN interactions between unique sites.
dists (dict) – NN, NNN, and NNNN interaction distances
ex_mat (DataFrame) – Invertible Heisenberg Hamiltonian for each graph.
ex_params (dict) – Exchange parameter values (meV/atom).
javg (float) – <J> exchange param (meV/atom).
igraph (StructureGraph) – Exchange interaction graph.
- class HeisenbergScreener(structures, energies, screen=False)[source]
Bases:
object
Clean and screen magnetic orderings.
Pre-processes magnetic orderings and energies for HeisenbergMapper. It prioritizes low-energy orderings with large and localized magnetic moments.
- Parameters:
structures (list) – Structure objects with magnetic moments.
energies (list) – Energies/atom of magnetic orderings.
screen (bool) – Try to screen out high energy and low-spin configurations.
pymatgen.analysis.magnetism.jahnteller module
JahnTeller distortion analysis.
- class JahnTellerAnalyzer[source]
Bases:
object
Will attempt to classify if structure may be Jahn-Teller active. Class currently uses datafile of hard-coded common Jahn-Teller active ions. If structure is annotated with magnetic moments, will estimate if structure may be high-spin or low-spin. Class aims for more false-positives than false-negatives.
Init for JahnTellerAnalyzer.
- get_analysis(structure: Structure, calculate_valences: bool = True, guesstimate_spin: bool = False, op_threshold: float = 0.1) dict [source]
Convenience method, uses get_analysis_and_structure method.
Obtain an analysis of a given structure and if it may be Jahn-Teller active or not. This is a heuristic, and may give false positives and false negatives (false positives are preferred).
- Parameters:
structure – input structure
calculate_valences – whether to attempt to calculate valences or not, structure should have oxidation states to perform analysis (Default value = True)
guesstimate_spin – whether to guesstimate spin state from magnetic moments or not, use with caution (Default value = False)
op_threshold – threshold for order parameter above which to consider site to match an octahedral or tetrahedral motif, since Jahn-Teller structures can often be quite distorted, this threshold is smaller than one might expect
- Returns:
analysis of structure, with key ‘strength’ which may be ‘none’, ‘strong’, ‘weak’, or ‘unknown’ (Default value = 0.1)
- get_analysis_and_structure(structure: Structure, calculate_valences: bool = True, guesstimate_spin: bool = False, op_threshold: float = 0.1) tuple[dict, Structure] [source]
Obtain an analysis of a given structure and if it may be Jahn-Teller active or not. This is a heuristic, and may give false positives and false negatives (false positives are preferred).
- Parameters:
structure – input structure
calculate_valences – whether to attempt to calculate valences or not, structure should have oxidation states to perform analysis (Default value = True)
guesstimate_spin – whether to guesstimate spin state from magnetic moments or not, use with caution (Default value = False)
op_threshold – threshold for order parameter above which to consider site to match an octahedral or tetrahedral motif, since Jahn-Teller structures can often be quite distorted, this threshold is smaller than one might expect
- Returns:
analysis of structure, with key ‘strength’ which may be ‘none’, ‘strong’, ‘weak’, or ‘unknown’ (Default value = 0.1) and decorated structure
- get_magnitude_of_effect_from_species(species: str | Species, spin_state: str, motif: str) str [source]
Get magnitude of Jahn-Teller effect from provided species, spin state and motif.
- Parameters:
species – e.g. Fe2+
spin_state – “high” or “low”
motif – “oct” or “tet”
- Returns:
“none”, “weak” or “strong”
- Return type:
str
- static get_magnitude_of_effect_from_spin_config(motif: str, spin_config: dict[str, float]) str [source]
Roughly, the magnitude of Jahn-Teller distortion will be: * in octahedral environments, strong if e_g orbitals unevenly occupied but weak if t_2g orbitals unevenly occupied * in tetrahedral environments always weaker.
- Parameters:
motif – “oct” or “tet”
spin_config – dict of ‘e’ (e_g) and ‘t’ (t2_g) with number of electrons in each state
- Returns:
“none”, “weak” or “strong”
- Return type:
str
- is_jahn_teller_active(structure: Structure, calculate_valences: bool = True, guesstimate_spin: bool = False, op_threshold: float = 0.1) bool [source]
Convenience method, uses get_analysis_and_structure method. Check if a given structure and if it may be Jahn-Teller active or not. This is a heuristic, and may give false positives and false negatives (false positives are preferred).
- Parameters:
structure – input structure
calculate_valences – whether to attempt to calculate valences or not, structure should have oxidation states to perform analysis (Default value = True)
guesstimate_spin – whether to guesstimate spin state from magnetic moments or not, use with caution (Default value = False)
op_threshold – threshold for order parameter above which to consider site to match an octahedral or tetrahedral motif, since Jahn-Teller structures can often be quite distorted, this threshold is smaller than one might expect
- Returns:
True if might be Jahn-Teller active, False if not
- Return type:
bool
- static mu_so(species: str | Species, motif: Literal['oct', 'tet'], spin_state: Literal['high', 'low']) float | None [source]
Calculate the spin-only magnetic moment for a given species. Only supports transition metals.
- Parameters:
species – Species
motif ("oct" | "tet") – Tetrahedron or octahedron crystal site coordination
spin_state ("low" | "high") – Whether the species is in a high or low spin state
- Returns:
- Spin-only magnetic moment in Bohr magnetons or None if
species crystal field not defined
- Return type:
float
- tag_structure(structure: Structure, calculate_valences: bool = True, guesstimate_spin: bool = False, op_threshold: float = 0.1) Structure [source]
Convenience method, uses get_analysis_and_structure method. Add a “possible_jt_active” site property on Structure.
- Parameters:
structure – input structure
calculate_valences – whether to attempt to calculate valences or not, structure should have oxidation states to perform analysis (Default value = True)
guesstimate_spin – whether to guesstimate spin state from magnetic moments or not, use with caution (Default value = False)
op_threshold – threshold for order parameter above which to consider site to match an octahedral or tetrahedral motif, since Jahn-Teller structures can often be quite distorted, this threshold is smaller than one might expect
- Returns:
Decorated Structure, will be in primitive setting.