Skip to content

Commit

Permalink
StructureData: Fix the pbc constraints of get_pymatgen_structure (#…
Browse files Browse the repository at this point in the history
…6281)

The pbc for getting a `pymatgen.core.structure.Structure` class via the
`StructureData.get_pymatgen_structure` were constrained to be 3D, i.e.
`(True, True, True)`. The core module of pymatgen actually allows for
any pbc condition, and this is now fixed by using the
`pymatgen.core.structure.Lattice` class to specify, along the cell
information, also the pbc.
  • Loading branch information
bastonero authored Feb 12, 2024
1 parent 5aae874 commit adcce4b
Show file tree
Hide file tree
Showing 2 changed files with 56 additions and 26 deletions.
31 changes: 17 additions & 14 deletions src/aiida/orm/nodes/data/structure.py
Original file line number Diff line number Diff line change
Expand Up @@ -1219,9 +1219,8 @@ def get_ase(self):
return self._get_object_ase()

def get_pymatgen(self, **kwargs):
"""Get pymatgen object. Returns Structure for structures with
periodic boundary conditions (in three dimensions) and Molecule
otherwise.
"""Get pymatgen object. Returns pymatgen Structure for structures with periodic boundary conditions
(in 1D, 2D, 3D) and Molecule otherwise.
:param add_spin: True to add the spins to the pymatgen structure.
Default is False (no spin added).
Expand All @@ -1237,7 +1236,7 @@ def get_pymatgen(self, **kwargs):
return self._get_object_pymatgen(**kwargs)

def get_pymatgen_structure(self, **kwargs):
"""Get the pymatgen Structure object.
"""Get the pymatgen Structure object with any PBC, provided the cell is not singular.
:param add_spin: True to add the spins to the pymatgen structure.
Default is False (no spin added).
Expand All @@ -1253,8 +1252,8 @@ def get_pymatgen_structure(self, **kwargs):
:return: a pymatgen Structure object corresponding to this
:py:class:`StructureData <aiida.orm.nodes.data.structure.StructureData>`
object.
:raise ValueError: if periodic boundary conditions do not hold
in at least one dimension of real space.
:raise ValueError: if the cell is singular, e.g. when it has not been set.
Use `get_pymatgen_molecule` instead, or set a proper cell.
"""
return self._get_object_pymatgen_structure(**kwargs)

Expand Down Expand Up @@ -1741,7 +1740,7 @@ def _get_object_pymatgen(self, **kwargs):
.. note:: Requires the pymatgen module (version >= 3.0.13, usage
of earlier versions may cause errors).
"""
if self.pbc == (True, True, True):
if any(self.pbc):
return self._get_object_pymatgen_structure(**kwargs)

return self._get_object_pymatgen_molecule(**kwargs)
Expand All @@ -1762,21 +1761,21 @@ def _get_object_pymatgen_structure(self, **kwargs):
:return: a pymatgen Structure object corresponding to this
:py:class:`StructureData <aiida.orm.nodes.data.structure.StructureData>`
object
:raise ValueError: if periodic boundary conditions does not hold
in at least one dimension of real space; if there are partial occupancies
together with spins (defined by kind names ending with '1' or '2').
:raise ValueError: if the cell is not set (i.e. is the default one);
if there are partial occupancies together with spins
(defined by kind names ending with '1' or '2').
.. note:: Requires the pymatgen module (version >= 3.0.13, usage
of earlier versions may cause errors)
"""
from pymatgen.core.lattice import Lattice
from pymatgen.core.structure import Structure

if self.pbc != (True, True, True):
raise ValueError('Periodic boundary conditions must apply in all three dimensions of real space')

species = []
additional_kwargs = {}

lattice = Lattice(matrix=self.cell, pbc=self.pbc)

if kwargs.pop('add_spin', False) and any(n.endswith('1') or n.endswith('2') for n in self.get_kind_names()):
# case when spins are defined -> no partial occupancy allowed
from pymatgen.core.periodic_table import Specie
Expand Down Expand Up @@ -1813,7 +1812,11 @@ def _get_object_pymatgen_structure(self, **kwargs):
raise ValueError(f'Unrecognized parameters passed to pymatgen converter: {kwargs.keys()}')

positions = [list(x.position) for x in self.sites]
return Structure(self.cell, species, positions, coords_are_cartesian=True, **additional_kwargs)

try:
return Structure(lattice, species, positions, coords_are_cartesian=True, **additional_kwargs)
except ValueError as err:
raise ValueError('Singular cell detected. Probably the cell was not set?') from err

def _get_object_pymatgen_molecule(self, **kwargs):
"""Converts
Expand Down
51 changes: 39 additions & 12 deletions tests/test_dataclasses.py
Original file line number Diff line number Diff line change
Expand Up @@ -2269,20 +2269,47 @@ class TestPymatgenFromStructureData:
"""

@skip_pymatgen
def test_1(self):
"""Tests the check of periodic boundary conditions."""
struct = StructureData()
struct.set_cell([[1, 0, 0], [0, 1, 2], [3, 4, 5]])
struct.pbc = [True, True, True]
struct.get_pymatgen_structure()
@pytest.mark.parametrize(
'pbc', ((True, True, True), (True, True, False), (True, False, False), (False, False, False))
)
def test_get_pymatgen_structure_pbc(self, pbc):
"""Tests the check of periodic boundary conditions when using the `get_pymatgen_structure` method."""
import numpy as np

struct.pbc = [True, True, False]
with pytest.raises(ValueError):
struct.get_pymatgen_structure()
cell = np.diag((1, 1, 1)).tolist()
symbols = ['Ba', 'Ba', 'Zr', 'Zr', 'O', 'O', 'O', 'O', 'O', 'O']
structure = StructureData(cell=cell, pbc=pbc)

for symbol in symbols:
structure.append_atom(name=symbol, symbols=[symbol], position=[0, 0, 0])

pymatgen = structure.get_pymatgen_structure()
assert pymatgen.lattice.pbc == pbc
assert pymatgen.lattice.matrix.tolist() == cell

@skip_pymatgen
def test_no_pbc(self):
"""Tests the `get_pymatgen*` methods for a 0D system, i.e. no periodic boundary conditions.
We expect `get_pymatgen` to return a `pymatgen.core.structure.Molecule` instance, and
`get_pymatgen_structure` to except since we did not set a cell, i.e. singular matrix is given.
"""
from pymatgen.core.structure import Molecule

symbol, pbc = 'Si', (False, False, False)
structure = StructureData(pbc=pbc)
structure.append_atom(name=symbol, symbols=[symbol], position=[0, 0, 0])

pymatgen = structure.get_pymatgen()
assert isinstance(pymatgen, Molecule)

match = 'Singular cell detected. Probably the cell was not set?'
with pytest.raises(ValueError, match=match):
structure.get_pymatgen_structure()

@skip_ase
@skip_pymatgen
def test_2(self):
def test_roundtrip_ase_aiida_pymatgen_structure(self):
"""Tests ASE -> StructureData -> pymatgen."""
import ase

Expand All @@ -2308,7 +2335,7 @@ def test_2(self):

@skip_ase
@skip_pymatgen
def test_3(self):
def test_roundtrip_ase_aiida_pymatgen_molecule(self):
"""Tests the conversion of StructureData to pymatgen's Molecule
(ASE -> StructureData -> pymatgen)
"""
Expand Down Expand Up @@ -2336,7 +2363,7 @@ def test_3(self):
]

@skip_pymatgen
def test_roundtrip(self):
def test_roundtrip_aiida_pymatgen_aiida(self):
"""Tests roundtrip StructureData -> pymatgen -> StructureData
(no spins)
"""
Expand Down

0 comments on commit adcce4b

Please sign in to comment.