Representation of molecules and crystals

The Structure class is designed to contain all necessary details of a molecular or crystalline structure. Moreover, it contains all methods that operate on a single structure, e.g. order parameters or the distance between atoms, and methods to manipulate the Structure object, e.g. deleting atoms or scaling the structrue’s unit cell. In case of the manipulation methods, a new Structure instance will be returned and the initial Structure object remains unchanged.

Initialization of the object

A structure is defined by its periodic boundary condition (pbc) and two lists with the length equal to the number of atoms:

  • elements contains the element symbol or number of each atom.

  • positions contains the cartesian coordinates (list or tuple of three numbers) of each atom.

The periodic boundary conditions can be given as single boolean or a list of 3 booleans for mixed periodic and non-periodic boundary conditions. If pbc is set to True for one or more directions a cell needs to be specified as well via a 3 times 3 nested list. If cell is given, the positions can also be passed on as scaled coordinates by setting the is_cartesian key to False or wrapped back into the unit cell by setting wrap to True.

The Structure object is created by passing the information at the initialization of the class:

[1]:
from aim2dat.strct import Structure

strct_molecule = Structure(
    label="molecule",
    elements=["N", "H", 1, 1],
    positions=[
        [0.000000, 0.000000, 0.000000],
        [0.000000, 0.000000, 1.008000],
        [0.950353, 0.000000, -0.336000],
        [-0.475176, -0.823029, -0.336000],
    ],
    pbc=False,
    site_attributes = {"oxidation_state": [-3, 1, 1, 1]}
)

strct_crystal = Structure(
    elements=["Ga", "As"],
    positions=[
        [0.0, 0.0, 0.0],
        [1.75, 1.75, 1.75],
    ],
    cell=[
        [0.0, 4.066, 4.0660001],
        [4.066, 0.0, 4.066],
        [4.066, 4.066, 0.0],
    ],
    is_cartesian=False,
    wrap=True,
    pbc=[True, True, True],
)

The following additional properties can be specified during the instance creation:

  • cell defines the unit cell for periodic structures.

  • is_cartesian specifies whether the given coordinates (positions) are cartesian or scaled by the cell vectors.

  • kinds can be assigned to distinguish atoms of the same element but with different properties (e.g. spin-up and spin-down in magnetic materials). It expects a list with the same length as elements and positions.

  • attributes stores a dictionary of different structural properties or can contain additional information about the structure.

  • site_attributes stores attributes containing information on the atom sites. Therefore, each value must have the same length as elements and positions.

  • site_attributes stores a dictionary likewise to attributes. However, each value needs to have the length of sites of the structure and individual values can be accessed via the iter_sites.

  • extras contains a dictionary with extensive calculation results, e.g. a partial radial distribution function.

  • label is a string and used as a key or identifier in the StructureCollection object.

An overview of the set properties is given by using the string-representation of the object:

[2]:
print(strct_crystal)
----------------------------------------------------------------------
-------------------------- Structure: None ---------------------------
----------------------------------------------------------------------

 Formula: GaAs
 PBC: [True True True]

                                 Cell
 Vectors: - [  0.0000   4.0660   4.0660]
          - [  4.0660   0.0000   4.0660]
          - [  4.0660   4.0660   0.0000]
 Lengths: [  5.7502   5.7502   5.7502]
 Angles: [ 60.0000  60.0000  60.0000]
 Volume: 134.4411

                                Sites
  - Ga None  [  0.0000   0.0000   0.0000] [  0.0000   0.0000   0.0000]
  - As None  [  6.0990   6.0990   6.0990] [  0.7500   0.7500   0.7500]
----------------------------------------------------------------------

Interface to other packages

A Structure object can also be directly created via the import-methods with prefix from_*. A list of all methods can be given via the import_methods property:

[3]:
Structure.list_import_methods()
[3]:
['from_file',
 'from_ase_atoms',
 'from_pymatgen_structure',
 'from_aiida_structuredata']

As an example we create an ase Atoms object (more details are given in the ASE documentation and use the corresponding import method to initiate the Structure object:

[4]:
from ase.spacegroup import crystal

a = 4.066 * 2.0
GaAs_conv = crystal(
    ("Ga", "As"),
    basis=((0.0, 0.0, 0.0), (0.75, 0.75, 0.75)),
    spacegroup=216,
    cellpar=[a, a, a, 90, 90, 90],
    primitive_cell=False,
)

strct_crystal_conv = Structure.from_ase_atoms(GaAs_conv, label="GaAs_crystal")
print(strct_crystal_conv)
----------------------------------------------------------------------
---------------------- Structure: GaAs_crystal -----------------------
----------------------------------------------------------------------

 Formula: Ga4As4
 PBC: [True True True]

                                 Cell
 Vectors: - [  8.1320   0.0000   0.0000]
          - [  0.0000   8.1320   0.0000]
          - [  0.0000   0.0000   8.1320]
 Lengths: [  8.1320   8.1320   8.1320]
 Angles: [ 90.0000  90.0000  90.0000]
 Volume: 537.7645

                                Sites
  - Ga None  [  0.0000   0.0000   0.0000] [  0.0000   0.0000   0.0000]
  - Ga None  [  0.0000   4.0660   4.0660] [  0.0000   0.5000   0.5000]
  - Ga None  [  4.0660   0.0000   4.0660] [  0.5000   0.0000   0.5000]
  - Ga None  [  4.0660   4.0660   0.0000] [  0.5000   0.5000   0.0000]
  - As None  [  6.0990   6.0990   6.0990] [  0.7500   0.7500   0.7500]
  - As None  [  2.0330   2.0330   6.0990] [  0.2500   0.2500   0.7500]
  - As None  [  2.0330   6.0990   2.0330] [  0.2500   0.7500   0.2500]
  - As None  [  6.0990   2.0330   2.0330] [  0.7500   0.2500   0.2500]
----------------------------------------------------------------------

In a similar way a Structure can be exported into the representing python objects of other python packages or written into a file (using the ASE package as backend). A list of all methods is given via the export_methods:

[5]:
Structure.list_export_methods()
[5]:
['to_dict',
 'to_file',
 'to_ase_atoms',
 'to_pymatgen_structure',
 'to_aiida_structuredata']

For example, a file is created using the to_file method:

[6]:
strct_crystal_conv.to_file("GaAs_crystal.xsf")

Manipulating structures

Another set of methods implemented in the Structure class act on the structure itself:

[7]:
strct_molecule.list_manipulation_methods()
[7]:
['delete_atoms', 'scale_unit_cell', 'substitute_elements']

The function returns a new Structure in case structural changes are induces otherwise the initial Structure object is returned. In general, kinds, attributes and site_attributes are transferred to the new Structure object (if not speficied otherwise). E.g. if “Ga” should be exchanged with “Al” we can simply call the substitute_elements function and obtain a new Structure object:

[8]:
strct_crystal_AsAl = strct_crystal_conv.substitute_elements(elements=["Ga", "Al"])
print(strct_crystal_AsAl)
----------------------------------------------------------------------
---------------------- Structure: GaAs_crystal -----------------------
----------------------------------------------------------------------

 Formula: Al4As4
 PBC: [True True True]

                                 Cell
 Vectors: - [  8.0987   0.0000   0.0000]
          - [  0.0000   8.0987   0.0000]
          - [  0.0000   0.0000   8.0987]
 Lengths: [  8.0987   8.0987   8.0987]
 Angles: [ 90.0000  90.0000  90.0000]
 Volume: 531.1797

                                Sites
  - Al None  [  0.0000   0.0000   0.0000] [  0.0000   0.0000   0.0000]
  - Al None  [  0.0000   4.0493   4.0493] [  0.0000   0.5000   0.5000]
  - Al None  [  4.0493   0.0000   4.0493] [  0.5000   0.0000   0.5000]
  - Al None  [  4.0493   4.0493   0.0000] [  0.5000   0.5000   0.0000]
  - As None  [  6.0740   6.0740   6.0740] [  0.7500   0.7500   0.7500]
  - As None  [  2.0247   2.0247   6.0740] [  0.2500   0.2500   0.7500]
  - As None  [  2.0247   6.0740   2.0247] [  0.2500   0.7500   0.2500]
  - As None  [  6.0740   2.0247   2.0247] [  0.7500   0.2500   0.2500]
----------------------------------------------------------------------

If the element is not part of the structure, no substitution can be conducted and the original object is returned

[9]:
strct_crystal_GaAs = strct_crystal_conv.substitute_elements(elements=["C", "Al"])
print(
    "Both structure objects have the same id and remain unchanged:",
    id(strct_crystal_GaAs) == id(strct_crystal_conv),
)
print(strct_crystal_GaAs)
Both structure objects have the same id and remain unchanged: True
----------------------------------------------------------------------
---------------------- Structure: GaAs_crystal -----------------------
----------------------------------------------------------------------

 Formula: Ga4As4
 PBC: [True True True]

                                 Cell
 Vectors: - [  8.1320   0.0000   0.0000]
          - [  0.0000   8.1320   0.0000]
          - [  0.0000   0.0000   8.1320]
 Lengths: [  8.1320   8.1320   8.1320]
 Angles: [ 90.0000  90.0000  90.0000]
 Volume: 537.7645

                                Sites
  - Ga None  [  0.0000   0.0000   0.0000] [  0.0000   0.0000   0.0000]
  - Ga None  [  0.0000   4.0660   4.0660] [  0.0000   0.5000   0.5000]
  - Ga None  [  4.0660   0.0000   4.0660] [  0.5000   0.0000   0.5000]
  - Ga None  [  4.0660   4.0660   0.0000] [  0.5000   0.5000   0.0000]
  - As None  [  6.0990   6.0990   6.0990] [  0.7500   0.7500   0.7500]
  - As None  [  2.0330   2.0330   6.0990] [  0.2500   0.2500   0.7500]
  - As None  [  2.0330   6.0990   2.0330] [  0.2500   0.7500   0.2500]
  - As None  [  6.0990   2.0330   2.0330] [  0.7500   0.2500   0.2500]
----------------------------------------------------------------------

External manipulation methods

In addition to the methods implemented in the Structure class less frequenetly used methods are implemented as functions in the ext_manipulation module. To use these methods the Structure object is passed as the first argument.

It is also possible to define new manipulation methods by using the decorator external_manipulation_method. The first argument needs to be named structure with type Structure and the last one change_label with type bool. If structural changes are introduced the method should return either a Structure object or a dictionary with the keyword arguments to initialize a Structure object as well as the label-suffix as string object. Otherwise it should return Ǹone.

The function add_atom shown below gives an example for an external manipulation method. In case the element does not occur in the structure the original structure is returned. If an atom of the same element is already present in the structure the method adds an atom with the element and position defined by element and position, respectively:

[10]:
from aim2dat.strct.ext_manipulation.decorator import external_manipulation_method
from aim2dat.utils.element_properties import get_element_symbol

@external_manipulation_method
def add_atom(structure, element, position, change_label=False):
    element = get_element_symbol(element)
    if element in structure.elements:
        new_structure = structure.to_dict()
        new_structure["elements"] = list(new_structure["elements"]) + [element]
        new_structure["positions"] = list(new_structure["positions"]) + [position]
        if new_structure["kinds"] is not None:
            new_structure["kinds"] = list(new_structure["kinds"]) + [None]
        for site_attr, val in new_structure["site_attributes"].items():
            new_structure["site_attributes"][site_attr] = list(val) + [None]
        return new_structure, "_add"

new_structure = add_atom(strct_molecule, "O", [3.0, 3.0, 3.0], change_label=True)
print(new_structure)
new_structure = add_atom(strct_molecule, "H", [3.0, 3.0, 3.0], change_label=True)
print(new_structure)
----------------------------------------------------------------------
------------------------ Structure: molecule -------------------------
----------------------------------------------------------------------

 Formula: NH3
 PBC: [False False False]

                                Sites
  - N  None  [  0.0000   0.0000   0.0000]
  - H  None  [  0.0000   0.0000   1.0080]
  - H  None  [  0.9504   0.0000  -0.3360]
  - H  None  [ -0.4752  -0.8230  -0.3360]
----------------------------------------------------------------------
----------------------------------------------------------------------
---------------------- Structure: molecule_add -----------------------
----------------------------------------------------------------------

 Formula: NH4
 PBC: [False False False]

                                Sites
  - N  None  [  0.0000   0.0000   0.0000]
  - H  None  [  0.0000   0.0000   1.0080]
  - H  None  [  0.9504   0.0000  -0.3360]
  - H  None  [ -0.4752  -0.8230  -0.3360]
  - H  None  [  3.0000   3.0000   3.0000]
----------------------------------------------------------------------