Skip to content

Commit 190c3e7

Browse files
authored
Fix default transformation_kwargs in analysis.magnetism.analyze.MagneticStructureEnumerator (#4275)
* minor type and comment cleanup * fix default kwargs dict * add unit test * skip if no enumlib * refactor to one liner
1 parent 7198b86 commit 190c3e7

File tree

2 files changed

+88
-65
lines changed

2 files changed

+88
-65
lines changed

src/pymatgen/analysis/magnetism/analyzer.py

Lines changed: 69 additions & 59 deletions
Original file line numberDiff line numberDiff line change
@@ -26,9 +26,10 @@
2626
from pymatgen.util.due import Doi, due
2727

2828
if TYPE_CHECKING:
29-
from typing import Any
29+
from collections.abc import Sequence
30+
from typing import Any, ClassVar
3031

31-
from numpy.typing import ArrayLike
32+
from numpy.typing import ArrayLike, NDArray
3233

3334
__author__ = "Matthew Horton"
3435
__copyright__ = "Copyright 2017, The Materials Project"
@@ -88,7 +89,7 @@ def __init__(
8889
threshold: float = 0,
8990
threshold_nonmag: float = 0.1,
9091
threshold_ordering: float = 1e-8,
91-
):
92+
) -> None:
9293
"""
9394
If magnetic moments are not defined, moments will be
9495
taken either from default_magmoms.yaml (similar to the
@@ -148,7 +149,7 @@ def __init__(
148149

149150
structure = structure.copy()
150151

151-
# check for disorder
152+
# Check for disorder structure
152153
if not structure.is_ordered:
153154
raise NotImplementedError(
154155
f"{type(self).__name__} not implemented for disordered structures, make ordered approximation first."
@@ -285,8 +286,28 @@ def __init__(
285286
self.structure = structure
286287
self.threshold_ordering = threshold_ordering
287288

289+
def __str__(self) -> str:
290+
"""Sorts a Structure (by fractional coordinate), and
291+
prints sites with magnetic information. This is
292+
useful over Structure.__str__ because sites are in
293+
a consistent order, which makes visual comparison between
294+
two identical Structures with different magnetic orderings
295+
easier.
296+
"""
297+
frac_coords = self.structure.frac_coords
298+
sorted_indices = np.lexsort((frac_coords[:, 2], frac_coords[:, 1], frac_coords[:, 0]))
299+
struct = Structure.from_sites([self.structure[idx] for idx in sorted_indices])
300+
301+
# adapted from Structure.__repr__
302+
outs = ["Structure Summary", repr(struct.lattice)]
303+
outs.append("Magmoms Sites")
304+
for site in struct:
305+
prefix = f"{site.properties['magmom']:+.2f} " if site.properties["magmom"] != 0 else " "
306+
outs.append(prefix + repr(site))
307+
return "\n".join(outs)
308+
288309
@staticmethod
289-
def _round_magmoms(magmoms: ArrayLike, round_magmoms_mode: float) -> np.ndarray:
310+
def _round_magmoms(magmoms: ArrayLike, round_magmoms_mode: float) -> NDArray:
290311
"""If round_magmoms_mode is an integer, simply round to that number
291312
of decimal places, else if set to a float will try and round
292313
intelligently by grouping magmoms.
@@ -334,7 +355,10 @@ def get_structure_with_spin(self) -> Structure:
334355

335356
return structure
336357

337-
def get_structure_with_only_magnetic_atoms(self, make_primitive: bool = True) -> Structure:
358+
def get_structure_with_only_magnetic_atoms(
359+
self,
360+
make_primitive: bool = True,
361+
) -> Structure:
338362
"""Get a Structure with only magnetic atoms present.
339363
340364
Args:
@@ -393,12 +417,12 @@ def get_ferromagnetic_structure(self, make_primitive: bool = True) -> Structure:
393417

394418
@property
395419
def is_magnetic(self) -> bool:
396-
"""Convenience property, returns True if any non-zero magmoms present."""
420+
"""True if any non-zero magmoms present."""
397421
return any(map(abs, self.structure.site_properties["magmom"]))
398422

399423
@property
400-
def magmoms(self) -> np.ndarray:
401-
"""Convenience property, returns magmoms as a numpy array."""
424+
def magmoms(self) -> NDArray:
425+
"""magmoms as a NumPy array."""
402426
return np.array(self.structure.site_properties["magmom"])
403427

404428
@property
@@ -448,11 +472,15 @@ def number_of_magnetic_sites(self) -> int:
448472
"""Number of magnetic sites present in structure."""
449473
return int(np.sum([abs(m) > 0 for m in self.magmoms]))
450474

451-
def number_of_unique_magnetic_sites(self, symprec: float = 1e-3, angle_tolerance: float = 5) -> int:
475+
def number_of_unique_magnetic_sites(
476+
self,
477+
symprec: float = 1e-3,
478+
angle_tolerance: float = 5,
479+
) -> int:
452480
"""
453481
Args:
454-
symprec: same as in SpacegroupAnalyzer (Default value = 1e-3)
455-
angle_tolerance: same as in SpacegroupAnalyzer (Default value = 5).
482+
symprec (float): same as in SpacegroupAnalyzer (Default value = 1e-3)
483+
angle_tolerance (float): same as in SpacegroupAnalyzer (Default value = 5).
456484
457485
Returns:
458486
int: Number of symmetrically-distinct magnetic sites present in structure.
@@ -509,7 +537,11 @@ def ordering(self) -> Ordering:
509537
return Ordering.AFM
510538
return Ordering.NM
511539

512-
def get_exchange_group_info(self, symprec: float = 1e-2, angle_tolerance: float = 5) -> tuple[str, int]:
540+
def get_exchange_group_info(
541+
self,
542+
symprec: float = 1e-2,
543+
angle_tolerance: float = 5,
544+
) -> tuple[str, int]:
513545
"""Get the information on the symmetry of the Hamiltonian
514546
describing the exchange energy of the system, taking into
515547
account relative direction of magnetic moments but not their
@@ -521,11 +553,11 @@ def get_exchange_group_info(self, symprec: float = 1e-2, angle_tolerance: float
521553
orderings within pymatgen.
522554
523555
Args:
524-
symprec: same as SpacegroupAnalyzer (Default value = 1e-2)
525-
angle_tolerance: same as SpacegroupAnalyzer (Default value = 5)
556+
symprec (float): same as SpacegroupAnalyzer (Default value = 1e-2)
557+
angle_tolerance (float): same as SpacegroupAnalyzer (Default value = 5)
526558
527559
Returns:
528-
spacegroup_symbol, international_number
560+
tuple[str, int]: spacegroup_symbol, international_number
529561
"""
530562
structure = self.get_structure_with_spin()
531563

@@ -562,34 +594,14 @@ def matches_ordering(self, other: Structure) -> bool:
562594

563595
return cmag_analyzer.matches(b_positive) or cmag_analyzer.matches(analyzer)
564596

565-
def __str__(self):
566-
"""Sorts a Structure (by fractional coordinate), and
567-
prints sites with magnetic information. This is
568-
useful over Structure.__str__ because sites are in
569-
a consistent order, which makes visual comparison between
570-
two identical Structures with different magnetic orderings
571-
easier.
572-
"""
573-
frac_coords = self.structure.frac_coords
574-
sorted_indices = np.lexsort((frac_coords[:, 2], frac_coords[:, 1], frac_coords[:, 0]))
575-
struct = Structure.from_sites([self.structure[idx] for idx in sorted_indices])
576-
577-
# adapted from Structure.__repr__
578-
outs = ["Structure Summary", repr(struct.lattice)]
579-
outs.append("Magmoms Sites")
580-
for site in struct:
581-
prefix = f"{site.properties['magmom']:+.2f} " if site.properties["magmom"] != 0 else " "
582-
outs.append(prefix + repr(site))
583-
return "\n".join(outs)
584-
585597

586598
class MagneticStructureEnumerator:
587599
"""Combines MagneticStructureAnalyzer and MagOrderingTransformation to
588600
automatically generate a set of transformations for a given structure
589601
and produce a list of plausible magnetic orderings.
590602
"""
591603

592-
available_strategies = (
604+
available_strategies: ClassVar[tuple[str, ...]] = (
593605
"ferromagnetic",
594606
"antiferromagnetic",
595607
"ferrimagnetic_by_motif",
@@ -602,14 +614,14 @@ def __init__(
602614
self,
603615
structure: Structure,
604616
default_magmoms: dict[str, float] | None = None,
605-
strategies: list[str] | tuple[str, ...] = (
617+
strategies: Sequence[str] = (
606618
"ferromagnetic",
607619
"antiferromagnetic",
608620
),
609621
automatic: bool = True,
610622
truncate_by_symmetry: bool = True,
611623
transformation_kwargs: dict | None = None,
612-
):
624+
) -> None:
613625
"""Generate different collinear magnetic orderings for a given input structure.
614626
615627
If the input structure has magnetic moments defined, it
@@ -618,18 +630,18 @@ def __init__(
618630
(this can be changed using default_magmoms kwarg).
619631
620632
Args:
621-
structure: input structure
622-
default_magmoms: (optional, defaults provided) dict of
623-
magnetic elements to their initial magnetic moments in μB, generally
624-
these are chosen to be high-spin since they can relax to a low-spin
625-
configuration during a DFT electronic configuration
626-
strategies: different ordering strategies to use, choose from:
633+
structure (Structure): input structure
634+
default_magmoms (dict[str, float], optional): magnetic elements to their
635+
initial magnetic moments in μB, generally these are chosen to be
636+
high-spin since they can relax to a low-spin configuration during
637+
a DFT electronic configuration
638+
strategies (Sequence[str]): ordering strategies to use, choose from:
627639
ferromagnetic, antiferromagnetic, antiferromagnetic_by_motif,
628640
ferrimagnetic_by_motif and ferrimagnetic_by_species (here, "motif",
629641
means to use a different ordering parameter for symmetry inequivalent
630642
sites)
631-
automatic: if True, will automatically choose sensible strategies
632-
truncate_by_symmetry: if True, will remove very unsymmetrical
643+
automatic (bool): if True, will automatically choose sensible strategies
644+
truncate_by_symmetry (bool): if True, will remove very unsymmetrical
633645
orderings that are likely physically implausible
634646
transformation_kwargs: keyword arguments to pass to
635647
MagOrderingTransformation, to change automatic cell size limits, etc.
@@ -654,10 +666,7 @@ def __init__(
654666
self.max_unique_sites = 8
655667

656668
# kwargs to pass to transformation (ultimately to enumlib)
657-
default_transformation_kwargs = {"check_ordered_symmetry": False, "timeout": 5}
658-
transformation_kwargs = transformation_kwargs or {}
659-
transformation_kwargs.update(default_transformation_kwargs)
660-
self.transformation_kwargs = transformation_kwargs
669+
self.transformation_kwargs = {"check_ordered_symmetry": False, "timeout": 5} | (transformation_kwargs or {})
661670

662671
# our magnetically ordered structures will be
663672
# stored here once generated and also store which
@@ -729,7 +738,10 @@ def _sanitize_input_structure(struct: Structure) -> Structure:
729738

730739
return struct
731740

732-
def _generate_transformations(self, structure: Structure) -> dict[str, MagOrderingTransformation]:
741+
def _generate_transformations(
742+
self,
743+
structure: Structure,
744+
) -> dict[str, MagOrderingTransformation]:
733745
"""The central problem with trying to enumerate magnetic orderings is
734746
that we have to enumerate orderings that might plausibly be magnetic
735747
ground states, while not enumerating orderings that are physically
@@ -938,15 +950,13 @@ def _generate_ordered_structures(
938950
and self.ordered_structures_origins instance variables.
939951
940952
Args:
941-
sanitized_input_structure: A sanitized input structure
942-
(_sanitize_input_structure)
953+
sanitized_input_structure (Structure): A sanitized input structure
943954
transformations: A dict of transformations (values) and name of
944-
enumeration strategy (key), the enumeration strategy name is just
945-
for record keeping
946-
955+
enumeration strategy (key), the enumeration strategy name is just
956+
for record keeping
947957
948958
Returns:
949-
list[Structures]
959+
tuple[list[Structure], list[str]]
950960
"""
951961
ordered_structures = self.ordered_structures
952962
ordered_structures_origins = self.ordered_structure_origins
@@ -1071,7 +1081,7 @@ def magnetic_deformation(structure_A: Structure, structure_B: Structure) -> Magn
10711081
ferromagnetic structures.
10721082
10731083
Adapted from Bocarsly et al. 2017,
1074-
doi: 10.1021/acs.chemmater.6b04729
1084+
DOI: 10.1021/acs.chemmater.6b04729
10751085
10761086
Args:
10771087
structure_A: Structure

tests/analysis/magnetism/test_analyzer.py

Lines changed: 19 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -19,9 +19,9 @@
1919

2020
TEST_DIR = f"{TEST_FILES_DIR}/analysis/magnetic_orderings"
2121

22-
enum_cmd = which("enum.x") or which("multienum.x")
23-
makestr_cmd = which("makestr.x") or which("makeStr.x") or which("makeStr.py")
24-
enumlib_present = enum_cmd and makestr_cmd
22+
ENUM_CMD = which("enum.x") or which("multienum.x")
23+
MAKESTR_CMD = which("makestr.x") or which("makeStr.x") or which("makeStr.py")
24+
ENUMLIB_PRESENT = ENUM_CMD and MAKESTR_CMD
2525

2626

2727
class TestCollinearMagneticStructureAnalyzer(TestCase):
@@ -259,10 +259,10 @@ def test_missing_spin(self):
259259
assert msa.structure.site_properties["magmom"] == [-5, 5, 0, 0]
260260

261261

262+
@pytest.mark.skipif(not ENUMLIB_PRESENT, reason="enumlib not present")
262263
class TestMagneticStructureEnumerator:
263-
@pytest.mark.skipif(not enumlib_present, reason="enumlib not present")
264264
def test_ordering_enumeration(self):
265-
# simple afm
265+
# simple AFM
266266
structure = Structure.from_file(f"{TEST_DIR}/LaMnO3.json")
267267
enumerator = MagneticStructureEnumerator(structure)
268268
assert enumerator.input_origin == "afm"
@@ -277,7 +277,7 @@ def test_ordering_enumeration(self):
277277
enumerator = MagneticStructureEnumerator(structure)
278278
assert enumerator.input_origin == "afm_by_Cr"
279279

280-
# afm requiring large cell size
280+
# AFM requiring large cell size
281281
# (enable for further development of workflow, too slow for CI)
282282

283283
# structure = Structure.from_file(f"{ref_dir}/CuO.json")
@@ -297,6 +297,19 @@ def test_ordering_enumeration(self):
297297
)
298298
assert enumerator.input_origin == "afm_by_motif_2a"
299299

300+
def test_default_transformation_kwargs(self):
301+
structure = Structure.from_file(f"{TEST_DIR}/LaMnO3.json")
302+
303+
# Make sure user input would not be overwritten by default values
304+
transformation_kwargs = {"timeout": 10, "check_ordered_symmetry": True}
305+
enumerator = MagneticStructureEnumerator(structure, transformation_kwargs=transformation_kwargs)
306+
assert enumerator.transformation_kwargs["timeout"] == 10
307+
assert enumerator.transformation_kwargs["check_ordered_symmetry"] is True
308+
309+
enumerator = MagneticStructureEnumerator(structure, transformation_kwargs=None)
310+
assert enumerator.transformation_kwargs["timeout"] == 5
311+
assert enumerator.transformation_kwargs["check_ordered_symmetry"] is False
312+
300313

301314
class TestMagneticDeformation:
302315
def test_magnetic_deformation(self):

0 commit comments

Comments
 (0)