Skip to content

Commit

Permalink
Make Supercell direct subclass of Crystal
Browse files Browse the repository at this point in the history
  • Loading branch information
LaurentRDC committed Jun 1, 2021
1 parent 4b0729f commit 5278f09
Show file tree
Hide file tree
Showing 4 changed files with 39 additions and 61 deletions.
5 changes: 5 additions & 0 deletions CHANGELOG.rst
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,11 @@ What's new

.. currentmodule:: crystals

1.3.1
-----

* The distinction between :class:`Supercell` and :class:`Crystal` no longer exists; :class:`Supercell` objects can be used everywhere a :class:`Crystal` is expected.

1.3.0
-----

Expand Down
85 changes: 30 additions & 55 deletions crystals/crystal.py
Original file line number Diff line number Diff line change
Expand Up @@ -96,7 +96,7 @@ def __hash__(self):

@property
def unitcell(self):
""" Generator of atoms forming the crystal unit cell. """
"""Generator of atoms forming the crystal unit cell."""
return super().__iter__()

@lru_cache(
Expand Down Expand Up @@ -284,7 +284,7 @@ def from_ase(cls, atoms, **kwargs):
)

def _spglib_cell(self):
""" Returns an array in spglib's cell format. """
"""Returns an array in spglib's cell format."""
# To get symmetry information, we only give spglib the unit cell atoms
# This way, all spglib-related methods (like symmetry()) will act on the unit cell only.
# This distinction is important for Crystal subclasses, like Supercell.
Expand Down Expand Up @@ -390,7 +390,20 @@ def supercell(self, n1, n2, n3):
cell : AtomicStructure
Iterable of `crystals.Atom` objects following the supercell dimensions.
"""
return Supercell(crystal=self, dimensions=(n1, n2, n3))
multicell = list()
for atm in self:
for factors in product(range(n1), range(n2), range(n3)):
fractional_offset = np.asarray(factors)
newatm = Atom(
element=atm.element,
coords=atm.coords_fractional + fractional_offset,
displacement=atm.displacement,
magmom=atm.magmom,
occupancy=atm.occupancy,
)
multicell.append(newatm)

return Supercell(unitcell=multicell, lattice_vectors=self.lattice_vectors)

def symmetry(self, symprec=1e-2, angle_tolerance=-1.0):
"""
Expand Down Expand Up @@ -500,7 +513,7 @@ def symmetry_operations(self, symprec=1e-2):
dataset = get_symmetry(cell=self._spglib_cell(), symprec=symprec)

def _to_affine(r, t):
""" Convert rotation and translation into single 4x4 affine transformation """
"""Convert rotation and translation into single 4x4 affine transformation"""
m = np.eye(4)
m[:3, :3] = r
m[:3, -1] = t
Expand Down Expand Up @@ -562,42 +575,42 @@ def pack(r, t):

@property
def international_symbol(self):
""" International Tables of Crystallography space-group short symbol. """
"""International Tables of Crystallography space-group short symbol."""
return self.symmetry()["international_symbol"]

@property
def international_full(self):
""" International Tables of Crystallography space-group full symbol. """
"""International Tables of Crystallography space-group full symbol."""
return self.symmetry()["international_full"]

@property
def hall_symbol(self):
""" Hall symbol. """
"""Hall symbol."""
return self.symmetry()["hall_symbol"]

@property
def hm_symbol(self):
""" Hermann-Mauguin symbol. """
"""Hermann-Mauguin symbol."""
return self.symmetry()["hm_symbol"]

@property
def pointgroup(self):
""" International Tables of Crystallography point-group. """
"""International Tables of Crystallography point-group."""
return self.symmetry()["pointgroup"]

@property
def international_number(self):
""" International Tables of Crystallography space-group number (between 1 and 230). """
"""International Tables of Crystallography space-group number (between 1 and 230)."""
return self.symmetry()["international_number"]

@property
def hall_number(self):
""" Hall number (between 1 and 531). """
"""Hall number (between 1 and 531)."""
return self.symmetry()["hall_number"]

@property
def centering(self):
""" Centering type of this crystals. """
"""Centering type of this crystals."""
return self.symmetry()["centering"]

def indexed_by(self, lattice):
Expand Down Expand Up @@ -630,11 +643,11 @@ def indexed_by(self, lattice):
)

def __str__(self):
""" String representation of this instance. Atoms may be omitted. """
"""String representation of this instance. Atoms may be omitted."""
return self._to_string(natoms=10)

def __repr__(self):
""" Verbose string representation of this instance. """
"""Verbose string representation of this instance."""
return self._to_string(natoms=len(self))

def _to_string(self, natoms):
Expand Down Expand Up @@ -738,55 +751,17 @@ def to_ase(self, **kwargs):
return ase_atoms(self)


class Supercell(AtomicStructure):
class Supercell(Crystal):
"""
The :class:`Supercell` class is a set-like container that represents a
supercell of crystalline structures.
It is recommended that you do not instantiate a :class:`Supercell` by hand, but rather
It is strongly recommended that you do not instantiate a :class:`Supercell` by hand, but rather
create a :class:`Crystal` object and use the :meth:`Crystal.supercell` method.
To iterate over all atoms in the supercell, use this object as an iterable.
To recover the underlying crystal, use the :attr:`Supercell.crystal` attribute.
Parameters
----------
crystal : Crystal
Crystal object from which the supercell is assembled.
dimensions : 3-tuple of ints
Number of cell repeats along the ``a1``, ``a2``, and ``a3`` directions. For example,
``(1, 1, 1)`` represents the trivial supercell.
"""

def __init__(self, crystal, dimensions, **kwargs):
self.crystal = crystal
self.dimensions = dimensions

n1, n2, n3 = self.dimensions

atoms = list()
for atm in crystal:
for factors in product(range(n1), range(n2), range(n3)):
fractional_offset = np.asarray(factors)
newatm = Atom(
element=atm.element,
coords=atm.coords_fractional + fractional_offset,
lattice=Lattice(self.crystal.lattice_vectors),
displacement=atm.displacement,
magmom=atm.magmom,
occupancy=atm.occupancy,
)
atoms.append(newatm)

super().__init__(atoms=atoms)

def __repr__(self):
n1, n2, n3 = self.dimensions
preamble = f"< Supercell object with dimensions ({n1} x {n2} x {n3}) and the following unit cell:\n"

lines = repr(self.crystal).splitlines(keepends=True)
lines[0] = preamble
return "".join(lines)
pass


@unique
Expand Down
5 changes: 3 additions & 2 deletions crystals/tests/test_crystal.py
Original file line number Diff line number Diff line change
Expand Up @@ -194,9 +194,10 @@ def test_from_cod_new_dir():
@pytest.mark.parametrize("name", islice(Crystal.builtins, 20))
def test_supercell_constructors(name):
""" Test Supercell constructors for varyous 'builtin' structures """
s = Crystal.from_database(name).supercell(2, 2, 2)
c = Crystal.from_database(name)
s = c.supercell(2, 2, 2)

assert len(s) == 8 * len(s.crystal)
assert len(s) == 8 * len(c)


def test_indexed_by_trivial_reindexing():
Expand Down
5 changes: 1 addition & 4 deletions docs/guides/crystals.rst
Original file line number Diff line number Diff line change
Expand Up @@ -562,10 +562,7 @@ Let's imagine we want to create a 2x2x2 supercell of graphite:
Supercells are generated by copying unit cells along crystallographic axes. For example, a 2x3x5 supercell will be extended by 2x along
the ``a1`` lattice vector, extended by 3x along the ``a2`` lattice vector, and extended by 5x along the ``a3`` lattice vector.

Users can recover the :class:`Crystal` object underlying a supercell via the :attr:`Supercell.crystal` attribute:

>>> spcl = Crystal.from_database('C').supercell(2,3,4)
>>> graphite = spcl.crystal
:class:`Supercell` objects can be used anywhere a :class:`Crystal` object can be used.

Compatibility with ASE
----------------------
Expand Down

0 comments on commit 5278f09

Please sign in to comment.