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(
    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,
)

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.

  • 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.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.export_methods
[5]:
['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")

Structural analysis methods

The Structure class also offers a few methods to analyse the structure. A list is given via the analysis_methods property:

[7]:
strct_molecule.analysis_methods
[7]:
['determine_point_group',
 'determine_space_group',
 'calculate_distance',
 'calculate_angle',
 'calculate_dihedral_angle',
 'calculate_voronoi_tessellation',
 'calculate_coordination',
 'calculate_ffingerprint']

By calling the function the analysis is performed and the results returned. As an example we determine the space group of the crystal using the spglib package as backend:

[8]:
strct_crystal.determine_space_group()
[8]:
{'space_group': {'int_symbol': 'F-43m',
  'sg_number': 216,
  'hall_number': 512,
  'point_group_symbol': '-43m',
  'schoenflies_symbol': 'Td^2',
  'centrosymmetric': False},
 'equivalent_sites': [{'wyckoff': 'a', 'symmetry': '-43m', 'sites': [0]},
  {'wyckoff': 'c', 'symmetry': '-43m', 'sites': [1]}]}

The analysis methods store two different kind of output data:

  • If the property store_calculated_properties is set to True a full set of the data is stored in the extras property and returned by the function.

  • A representative subset of the data is stored in the attributes property.

[9]:
print("Attribute:", strct_crystal.attributes["space_group"])
print("Extra:", strct_crystal.extras["space_group"])
Attribute: 216
Extra: {'space_group': {'int_symbol': 'F-43m', 'sg_number': 216, 'hall_number': 512, 'point_group_symbol': '-43m', 'schoenflies_symbol': 'Td^2', 'centrosymmetric': False}, 'equivalent_sites': [{'wyckoff': 'a', 'symmetry': '-43m', 'sites': [0]}, {'wyckoff': 'c', 'symmetry': '-43m', 'sites': [1]}]}

In case the same analysis function of the instance with the same parameters is called again the stored data in the extras property is returned instead of performing the same analysis again. This feature is very handy but also demands sufficient memory capacities, it can be switched on and off via the store_calculated_properties property.

Manipulating structures

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

[10]:
strct_molecule.manipulation_methods
[10]:
['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. E.g. if “Ga” should be exchanged with “Al” we can simply call the substitute_elements function and obtain a new Structure object:

[11]:
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

[12]:
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 analysis and manipulation methods

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

[13]:
from aim2dat.strct.ext_analysis import calculate_prdf

prdf = calculate_prdf(strct_crystal_GaAs, r_max=7.5)

The calculated properties are still stored within the Structure object if the store_calculated_properties is set to True:

[14]:
strct_crystal_GaAs.extras.keys()
[14]:
dict_keys(['prdf'])

It is also possible to define individual _analysis_ or _manipulation_ methods by using the corresponding decorators external_analysis_method and external_manipulation_method, respectively.

For both, analysis and manipulation methods, the first argument needs to be the Structure object. The output for the analysis methods needs to be a tuple with a length of two consisting of calculated attributes and extras. The manipulation methods need to return a dictionary with the keyword arguments to initialize a Structure object if a manipulation is performed and otherwise Ǹone.

The function calculate_n_element presented below presents an example for an external analysis method. It determines how many atoms of a specific element are present in the Structure and returns it in extras.

[15]:
from aim2dat.strct.ext_analysis.decorator import external_analysis_method
from aim2dat.utils.element_properties import get_element_symbol


@external_analysis_method
def determine_n_element(structure, element="H"):
    element = get_element_symbol(element)
    if element in structure.chem_formula:
        return (None, structure.chem_formula[element])
    else:
        return (None, 0)


determine_n_element(strct_crystal_GaAs, "Ga")
strct_crystal_GaAs.extras["n_element"]
[15]:
4