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: pymatgen.core.structure.Structure, overwrite_magmom_mode: Union[pymatgen.analysis.magnetism.analyzer.OverwriteMagmomMode, str] = 'none', round_magmoms: bool = False, detect_valences: bool = False, make_primitive: bool = True, default_magmoms: dict = None, set_net_positive: bool = True, threshold: float = 0.0, threshold_nonmag: float = 0.1)[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 < theshold 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

get_exchange_group_info(symprec: float = 0.01, angle_tolerance: float = 5.0) → Tuple[str, int][source]

Returns 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.0)

Returns

spacegroup_symbol, international_number

get_ferromagnetic_structure(make_primitive: bool = True) → pymatgen.core.structure.Structure[source]

Returns 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) → pymatgen.core.structure.Structure[source]

Returns 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) → pymatgen.core.structure.Structure[source]

Returns 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() → pymatgen.core.structure.Structure[source]

Returns a Structure with species decorated with spin values instead of using magmom site properties.

property is_magnetic

Convenience property, returns True if any non-zero magmoms present.

property magmoms

Convenience property, returns magmoms as a numpy array.

property magnetic_species_and_magmoms

Returns 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: pymatgen.core.structure.Structure) → bool[source]

Compares the magnetic orderings of one structure with another.

Parameters

other – Structure to compare

Returns: True or False

property number_of_magnetic_sites

Number of magnetic sites present in structure.

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.

property ordering

Applies heuristics to return a magnetic ordering for a collinear magnetic structure. Result is not guaranteed for correctness.

Returns: Ordering Enum (‘FiM’ is used as the abbreviation for ferrimagnetic)

property types_of_magnetic_specie

Equivalent to Structure.types_of_specie but only returns magnetic species.

Returns: types of Specie as a list

class MagneticDeformation(type, deformation)

Bases: tuple

Create new instance of MagneticDeformation(type, deformation)

property deformation

Alias for field number 1

property type

Alias for field number 0

class MagneticStructureEnumerator(structure: pymatgen.core.structure.Structure, default_magmoms: Optional[Dict[str, float]] = None, strategies: Union[List[str], Tuple[str, ...]] = ('ferromagnetic', 'antiferromagnetic'), automatic: bool = True, truncate_by_symmetry: bool = True, transformation_kwargs: Optional[Dict] = 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.

This class will try generated 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.

available_strategies = ('ferromagnetic', 'antiferromagnetic', 'ferrimagnetic_by_motif', 'ferrimagnetic_by_species', 'antiferromagnetic_by_motif', 'nonmagnetic')
class Ordering[source]

Bases: enum.Enum

Enumeration defining possible magnetic orderings.

AFM = 'AFM'
FM = 'FM'
FiM = 'FiM'
NM = 'NM'
Unknown = 'Unknown'
class OverwriteMagmomMode[source]

Bases: enum.Enum

Enumeration defining different modes for analyzer.

none = 'none'
normalize = 'normalize'
replace_all = 'replace_all'
respect_sign = 'respect_sign'
respect_zero = 'respect_zeros'
magnetic_deformation(structure_A: pymatgen.core.structure.Structure, structure_B: pymatgen.core.structure.Structure) → pymatgen.analysis.magnetism.analyzer.MagneticDeformation[source]

Calculates ‘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: Magnetic deformation