aiida_crystal17.symmetry package

Submodules

aiida_crystal17.symmetry.symmetry module

A module for computing the symmetry of an AiiDA StructureData object.

When computing symmetry, atomic sites with the same Kind are treated as symmetrically equivalent (rather than just the atomic elements).

Currently only 3D structures are considered.

NB: this module is not specific to CRYSTAL, and may be move to a separate package at a later date

aiida_crystal17.symmetry.symmetry.affine_to_operation(affine_matrix)[source]

create a flattened symmetry operation, from a 4x4 affine transformation matrix

Parameters

affine_matrix (list) – 4x4 affine transformation

Returns

representing symmetry operation as a flattened list; (r00, r01, r02, r10, r11, r12, r20, r21, r22, t0, t1, t2)

Return type

list

aiida_crystal17.symmetry.symmetry.cartesian_to_frac(lattice, ccoords)[source]

convert cartesian coordinates to fractional coordinates

Parameters
  • lattice (list) – 3x3 array of lattice vectors

  • ccoords (list) – Nx3 array of cartesian coordinate

Returns

Nx3 array of fractional coordinate

Return type

list

aiida_crystal17.symmetry.symmetry.compute_symmetry_dataset(structure, symprec, angle_tolerance, use_kinds=True)[source]

compute the symmetry of a Structure, with periodic boundary conditions in all axes, using spglib.

When computing symmetry, atomic sites with the same Kind are treated as symmetrically equivalent (rather than just the atomic elements).

Parameters
  • structure (aiida.StructureData) –

  • symprec (float) – Symmetry search tolerance in the unit of length.

  • angle_tolerance (float or None) – Symmetry search tolerance in the unit of angle degrees. If the value is negative or None, an internally optimized routine is used to judge symmetry.

  • use_kinds (bool) – if True use kind names to define inequivalent sites, else use symbols

Returns

dataset – spglib symmetry dataset

Return type

dict

aiida_crystal17.symmetry.symmetry.compute_symmetry_dict(structure, symprec, angle_tolerance, use_kinds=True)[source]

compute the symmetry of a Structure, with periodic boundary conditions in all axes, using spglib

When computing symmetry, atomic sites with the same Kind are treated as symmetrically equivalent (rather than just the atomic elements).

Parameters
  • structure (aiida.StructureData) –

  • symprec (float) – Symmetry search tolerance in the unit of length.

  • angle_tolerance (float or None) – Symmetry search tolerance in the unit of angle degrees. If the value is negative or None, an internally optimized routine is used to judge symmetry.

  • use_kinds (bool) – if True use kind names to define inequivalent sites, else use symbols

Returns

data required to create an AiiDA SymmetryData object

Return type

dict

aiida_crystal17.symmetry.symmetry.convert_structure(structure, out_type)[source]

convert an AiiDA, ASE or dict object to another type

Parameters
aiida_crystal17.symmetry.symmetry.find_primitive(structure, symprec, angle_tolerance)[source]

compute the primitive cell for an AiiDA structure

When computing symmetry, atomic sites with the same Kind are treated as symmetrically equivalent (rather than just the atomic elements).

Parameters
  • structure (aiida.StructureData) –

  • symprec (float) – Symmetry search tolerance in the unit of length.

  • angle_tolerance (float or None) – Symmetry search tolerance in the unit of angle degrees. If the value is negative, an internally optimized routine is used to judge symmetry.

Returns

Return type

aiida.StructureData

aiida_crystal17.symmetry.symmetry.frac_to_cartesian(lattice, fcoords)[source]

convert fractional coordinates to cartesian coordinates

Parameters
  • lattice (list) – 3x3 array of lattice vectors

  • fcoords (list) – Nx3 array of fractional coordinate

Returns

Nx3 array of cartesian coordinate

Return type

list

aiida_crystal17.symmetry.symmetry.get_crystal_system_name(sg_number)[source]

Get the crystal system for the structure from the space group number

Parameters

sg_number (int) – the spacegroup number

Returns

crystal_system – Crystal system for structure

Return type

str

aiida_crystal17.symmetry.symmetry.get_hall_number_from_symmetry(operations, basis='fractional', lattice=None, symprec=1e-05)[source]

obtain the Hall number from the symmetry operations

Parameters
  • operations (list) – Nx12 flattened list of symmetry operations

  • basis (str) – “fractional” or “cartesian”

Returns

Return type

int

aiida_crystal17.symmetry.symmetry.get_lattice_type_name(sg_number)[source]

Get the lattice type for the structure from the space group

This is the same as crystal system name, with the exception of the trigonal -> hexagonal or rhombohedral

Parameters

sg_number (int) – the spacegroup number

Returns

lattice_type – Crystal lattice for the structure

Return type

str

aiida_crystal17.symmetry.symmetry.operation_cart_to_frac(lattice, rotation, translation)[source]

convert a single symmetry operation from cartesian to fractional

Parameters
  • lattice (list) – 3x3 matrix of lattice vectors (a, b, c)

  • rotation (list) – 3x3 rotation matrix

  • translation (list) – 3x1 translation vector

Returns

  • rotation (list) – 3x3 rotation matrix

  • translation (list) – 3x1 translation vector

aiida_crystal17.symmetry.symmetry.operation_frac_to_cart(lattice, rotation, translation)[source]

convert a single symmetry operation from fractional to cartesian

Parameters
  • lattice (list) – 3x3 matrix of lattice vectors (a, b, c)

  • rotation (list) – 3x3 rotation matrix

  • translation (list) – 3x1 translation vector

Returns

  • rotation (list) – 3x3 rotation matrix

  • translation (list) – 3x1 translation vector

aiida_crystal17.symmetry.symmetry.operation_to_affine(operation)[source]

create a 4x4 affine transformation matrix, from a flattened symmetry operation

Parameters

operation (list) – representing symmetry operation as a flattened list; (r00, r01, r02, r10, r11, r12, r20, r21, r22, t0, t1, t2)

Returns

4x4 array

Return type

list

aiida_crystal17.symmetry.symmetry.operations_cart_to_frac(operations, lattice)[source]

convert a list of cartesian symmetry operations to fractional

Parameters
  • operations (list) – Nx9 array, representing each symmetry operation as a flattened list; (r00, r01, r02, r10, r11, r12, r20, r21, r22, t0, t1, t2)

  • lattice (list) – 3x3 matrix (a, b, c)

Returns

Nx9 array of operations

Return type

list

aiida_crystal17.symmetry.symmetry.operations_frac_to_cart(operations, lattice)[source]

convert a list of fractional symmetry operations to cartesian

Parameters
  • operations (list) – Nx9 array, representing each symmetry operation as a flattened list; (r00, r01, r02, r10, r11, r12, r20, r21, r22, t0, t1, t2)

  • lattice (list) – 3x3 matrix (a, b, c)

Returns

Nx9 array of operations

Return type

list

aiida_crystal17.symmetry.symmetry.prepare_for_spglib(structure, use_kinds=True)[source]

prepare an AiiDa Structure for parsing to spglib, labelling sites with the same Kind as equivalent

Parameters
  • structure (aiida.StructureData) –

  • use_kinds (bool) – if True use kind names to define inequivalent sites, else use symbols

Returns

  • tuple (cell) – (lattice, fcoords, inequivalent)

  • dict (int2kind_map) – maps integer values in inequivalent list to AiiDa Kind objects

aiida_crystal17.symmetry.symmetry.print_structure(structure, max_srows=None, round_dp=4)[source]

print a formatted string, with information about a StructureData cell and sites

Parameters
  • structure (aiida.StructureData) –

  • max_srows (None or int) – limit the number of site lines returned

  • round_dp (int) – round numbers to n decimal places

aiida_crystal17.symmetry.symmetry.reset_kind_names(structure, kind_names)[source]

reset the kind names (per site) of a StructureData node

Parameters
Returns

a cloned node

Return type

aiida.StructureData

Raises

AssertionError – if the kind_names are not compatible with the current sites

aiida_crystal17.symmetry.symmetry.standardize_cell(structure, symprec, angle_tolerance, to_primitive=False, no_idealize=False)[source]

compute the standardised cell for an AiiDA structure

When computing symmetry, atomic sites with the same Kind are treated as symmetrically equivalent (rather than just the atomic elements).

Parameters
  • structure (aiida.StructureData) –

  • to_primitive (bool) – If True, the standardized primitive cell is created.

  • no_idealize (bool) – If True, it is disabled to idealize lengths and angles of basis vectors and positions of atoms according to crystal symmetry.

  • symprec (float) – Symmetry search tolerance in the unit of length.

  • angle_tolerance (float or None) – Symmetry search tolerance in the unit of angle degrees. If the value is negative or None, an internally optimized routine is used to judge symmetry.

Returns

Return type

aiida.StructureData

aiida_crystal17.symmetry.symmetry.structure_info(structure, max_srows=None, round_dp=4)[source]

get a formatted string, with information about a StructureData cell and sites

Parameters
  • structure (aiida.StructureData) –

  • max_srows (None or int) – limit the number of site lines returned

  • round_dp (int) – round numbers to n decimal places

Returns

Return type

str

aiida_crystal17.symmetry.symmetry.structure_to_dict(structure)[source]

create a dictionary of structure properties per atom

Parameters

structure (aiida.StructureData) – the input structure

Returns

dictionary containing; lattice, atomic_numbers, ccoords, pbc, kinds, equivalent

Return type

dict

Module contents