pyiron.atomistics.structure.atoms module

class pyiron.atomistics.structure.atoms.Atoms(symbols=None, positions=None, numbers=None, tags=None, momenta=None, masses=None, magmoms=None, charges=None, scaled_positions=None, cell=None, pbc=None, celldisp=None, constraint=None, calculator=None, info=None, indices=None, elements=None, dimension=None, species=None, **qwargs)[source]

Bases: object

The Atoms class represents all the information required to describe a structure at the atomic scale. This class is written in such a way that is compatible with the ASE atoms class. Some of the functions in this module is based on the corresponding implementation in the ASE package

Parameters
  • elements (list/numpy.ndarray) – List of strings containing the elements or a list of atomistics.structure.periodic_table.ChemicalElement instances

  • numbers (list/numpy.ndarray) – List of atomic numbers of elements

  • symbols (list/numpy.ndarray) – List of chemical symbols

  • positions (list/numpy.ndarray) – List of positions

  • scaled_positions (list/numpy.ndarray) – List of scaled positions (relative coordinates)

  • pbc (list/numpy.ndarray/boolean) – Tells if periodic boundary conditions should be applied on the three axes

  • cell (list/numpy.ndarray instance) – A 3x3 array representing the lattice vectors of the structure

Note: Only one of elements/symbols or numbers should be assigned during initialization

indices

A list of size N which gives the species index of the structure which has N atoms

Type

numpy.ndarray

add_high_symmetry_path(path)[source]

Adds a new path to the dictionary of pathes for band structure calculations.

Parameters

path (dict) – dictionary of lists of tuples with start and end point. E.G. {“my_path”: [(‘Gamma’, ‘X’), (‘X’, ‘Y’)]}

add_high_symmetry_points(new_points)[source]

Adds new points to the dict of existing high symmetry points.

Parameters

new_points (dict) – Points to add

add_tag(*args, **qwargs)[source]

Add tags to the atoms object.

Examples

For selective dynamics:

>>> self.add_tag(selective_dynamics=[False, False, False])
analyse_ovito_centro_symmetry(num_neighbors=12)[source]
analyse_ovito_cna_adaptive(mode='total')[source]

Use Ovito’s common neighbor analysis binding.

Parameters

mode ("total"/"numeric"/"str") – Controls the style and level of detail of the output. (Default is “total”, only return a summary of the values in the structure.)

Returns

(depends on mode)

analyse_ovito_voronoi_volume()[source]
analyse_phonopy_equivalent_atoms()[source]
analyse_pyscal_steinhardt_parameter(cutoff=3.5, n_clusters=2, q=[4, 6])[source]
append(atom)[source]

Append atom to end. Copied from ase

Parameters

atom

Returns:

apply_strain(epsilon, return_box=False)[source]
Parameters
  • epsilon (float/list/ndarray) – epsilon matrix. If a single number is set, the same strain is applied in each direction. If a 3-dim vector is set, it will be multiplied by a unit matrix.

  • return_box (bool) – whether to return a box. If set to True, only the returned box will have the desired strain and the original box will stay unchanged.

property cell

A size 3x3 array which gives the lattice vectors of the cell as [a1, a2, a3]

Type

numpy.ndarray

center(vacuum=None, axis=0, 1, 2)[source]

Center atoms in unit cell.

Adopted from ASE code (https://wiki.fysik.dtu.dk/ase/_modules/ase/atoms.html#Atoms.center)

Parameters
  • vacuum (float) – If specified adjust the amount of vacuum when centering. If vacuum=10.0 there will thus be 10 Angstrom of vacuum on each side.

  • axis (tuple/list) – List or turple of integers specifying the axis along which the atoms should be centered

center_coordinates_in_unit_cell(origin=0, eps=0.0001)[source]

Wrap atomic coordinates within the supercell as given by a1, a2., a3

Parameters
  • origin (float) – 0 to confine between 0 and 1, -0.5 to confine between -0.5 and 0.5

  • eps (float) – Tolerance to detect atoms at cell edges

Returns

Wrapped structure

Return type

pyiron.atomistics.structure.atoms.Atoms

close()[source]
cluster_analysis(id_list, neighbors=None, radius=None, return_cluster_sizes=False)[source]
Parameters
  • id_list

  • neighbors

  • radius

  • return_cluster_sizes

Returns:

convert_element(el, pse=None)[source]

Convert a string or an atom instance into a ChemicalElement instance

Parameters
  • el (str/atomistics.structure.atom.Atom) – String or atom instance from which the element should be generated

  • pse (atomistics.structure.periodictable.PeriodicTable) – PeriodicTable instance from which the element is generated (optional)

Returns

The required chemical element

Return type

atomistics.structure.periodictable.ChemicalElement

static convert_formula(elements)[source]
Parameters

elements

Returns:

copy()[source]

Returns a copy of the instance

Returns

A copy of the instance

Return type

pyiron.atomistics.structure.atoms.Atoms

create_line_mode_structure(with_time_reversal=True, recipe='hpkot', threshold=1e-07, symprec=1e-05, angle_tolerance=- 1.0)[source]

Uses ‘seekpath’ to create a new structure with high symmetry points and path for band structure calculations.

Parameters
  • with_time_reversal (bool) – if False, and the group has no inversion symmetry, additional lines are returned as described in the HPKOT paper.

  • recipe (str) – choose the reference publication that defines the special points and paths. Currently, only ‘hpkot’ is implemented.

  • threshold (float) – the threshold to use to verify if we are in and edge case (e.g., a tetragonal cell, but a==c). For instance, in the tI lattice, if abs(a-c) < threshold, a EdgeCaseWarning is issued. Note that depending on the bravais lattice, the meaning of the threshold is different (angle, length, …)

  • symprec (float) – the symmetry precision used internally by SPGLIB

  • angle_tolerance (float) – the angle_tolerance used internally by SPGLIB

Returns

new structure

Return type

pyiron.atomistics.structure.atoms.Atoms

property elements

A size N list of atomistics.structure.periodic_table.ChemicalElement instances according to the ordering of the atoms in the instance

Type

numpy.ndarray

extend(other)[source]

Extend atoms object by appending atoms from other. Copied from ase

Parameters

other

Returns:

find_neighbors_by_vector(vector, deviation=False, num_neighbors=96)[source]
Parameters
  • vector (list/np.ndarray) – vector by which positions are translated (and neighbors are searched)

  • deviation (bool) – whether to return distance between the expect positions and real positions

  • num_neighbors (int) – number of neighbors to take into account in get_neighbors

Returns

list of id’s for the specified translation

Return type

np.ndarray

Example

a_0 = 2.832 structure = pr.create_structure(‘Fe’, ‘bcc’, a_0) id_list = structure.find_neighbors_by_vector([0, 0, a_0]) # In this example, you get a list of neighbor atom id’s at z+=a_0 for each atom. # This is particularly powerful for SSA when the magnetic structure has to be translated # in each direction.

from_hdf(hdf, group_name='structure')[source]

Retrieve the object from a HDF5 file

Parameters
Returns

The retrieved atoms class

Return type

pyiron_atomistic.structure.atoms.Atoms

get_array(name, copy=True)[source]

Get an array. This function is for the purpose of compatibility with the ASE package

Parameters
  • name (str) – Name of the required array

  • copy (bool) – True if a copy of the array is to be returned

Returns

An array of a copy of the array

get_atomic_numbers()[source]

Returns the atomic numbers of all the atoms in the structure

Returns

A list of atomic numbers

Return type

numpy.ndarray

get_bonds(radius=None, max_shells=None, prec=0.1, num_neighbors=20)[source]
Parameters
  • radius

  • max_shells

  • prec – minimum distance between any two clusters (if smaller considered to be single cluster)

  • num_neighbors

Returns:

get_boundary_region(dist)[source]

get all atoms in the boundary around the supercell which have a distance to the supercell boundary of less than dist

Parameters

dist

Returns:

get_calculator()[source]

Get currently attached calculator object.

get_cell(complete=False)[source]

Get the three unit cell vectors as a 3x3 ndarray.

get_center_of_mass()[source]
Returns

center of mass in A

Return type

com (float)

get_chemical_elements()[source]

Returns the list of chemical element instances

Returns

A list of chemical element instances

Return type

numpy.ndarray

get_chemical_formula()[source]

Returns the chemical formula of structure

Returns

The chemical formula as a string

Return type

str

get_chemical_indices()[source]

Returns the list of chemical indices as ordered in self.species

Returns

A list of chemical indices

Return type

numpy.ndarray

get_chemical_symbols()[source]

Returns the chemical symbols for all the atoms in the structure

Returns

A list of chemical symbols

Return type

numpy.ndarray

get_constraint()[source]
get_density()[source]

Returns the density in g/cm^3

Returns

Density of the structure

Return type

float

get_distance(a0, a1, mic=True, vector=False)[source]

Return distance between two atoms.

Use mic=True to use the Minimum Image Convention. vector=True gives the distance vector (from a0 to a1).

Parameters
  • a0 – position or atom ID

  • a1 – position or atom ID

  • mic – minimum image convention (True if periodic boundary conditions should be considered)

  • vector – True, if instead of distnce the vector connecting the two positions should be returned

Returns: distance or vectors in length unit

get_distance_matrix(mic=True, vector=False)[source]

Return distances between all atoms in a matrix. cf. get_distance

get_distances(a0=None, a1=None, mic=True, vector=False)[source]

Return distance matrix of every position in p1 with every position in p2

Parameters
  • a0 (numpy.ndarray/list) – Nx3 array of positions

  • a1 (numpy.ndarray/list) – Nx3 array of positions

  • mic (bool) – minimum image convention

  • vector (bool) – return vectors instead of distances

Returns

numpy.ndarray NxN if vector=False and NxNx3 if vector=True

if a1 is not set, it is assumed that distances between all positions in a0 are desired. a1 will be set to a0 in this case. if both a0 and a1 are not set, the distances between all atoms in the box are returned

Use mic to use the minimum image convention.

Learn more about get_distances from the ase website: https://wiki.fysik.dtu.dk/ase/ase/geometry.html#ase.geometry.get_distances

get_equivalent_points(points, use_magmoms=False, use_elements=True, symprec=1e-05, angle_tolerance=- 1.0)[source]
Parameters
  • points (list/ndarray) – 3d vector

  • use_magmoms (bool) – cf. get_symmetry()

  • use_elements (bool) – cf. get_symmetry()

  • symprec (float) – cf. get_symmetry()

  • angle_tolerance (float) – cf. get_symmetry()

Returns

array of equivalent points with respect to box symmetries

Return type

(ndarray)

get_equivalent_voronoi_vertices(return_box=False, minimum_dist=0.1, symprec=1e-05, angle_tolerance=- 1.0)[source]

This function gives the positions of spatially equivalent Voronoi vertices in lists, which most likely represent interstitial points or vacancies (along with other high symmetry points) Each list item contains an array of positions which are spacially equivalent. This function does not work if there are Hs atoms in the box

Parameters
  • return_box – True, if the box containing atoms on the positions of Voronoi vertices should be returned (which are represented by Hs atoms)

  • minimum_dist – Minimum distance between two Voronoi vertices to be considered as one

Returns: List of numpy array positions of spacially equivalent Voronoi vertices

get_high_symmetry_path()[source]

Path used for band structure calculations

Returns

dict of pathes with start and end points.

Return type

dict

get_high_symmetry_points()[source]

dictionary of high-symmetry points defined for this specific structure.

Returns

high_symmetry_points

Return type

dict

get_initial_magnetic_moments()[source]

Get array of initial magnetic moments.

Returns

numpy.array()

get_ir_reciprocal_mesh(mesh, is_shift=array([0, 0, 0], dtype=int32), is_time_reversal=True, symprec=1e-05)[source]
Parameters
  • mesh

  • is_shift

  • is_time_reversal

  • symprec

Returns:

get_majority_species(return_count=False)[source]

This function returns the majority species and their number in the box

Returns

number of atoms of the majority species, chemical symbol and chemical index

get_masses()[source]

Returns:

get_masses_dof()[source]

Returns:

get_neighborhood(position, num_neighbors=12, t_vec=True, include_boundary=True, tolerance=2, id_list=None, cutoff=None, cutoff_radius=None)[source]
Parameters
  • position – position in a box whose neighborhood information is analysed

  • num_neighbors

  • t_vec (bool) – True: compute distance vectors (pbc are automatically taken into account)

  • include_boundary (bool) – True: search for neighbors assuming periodic boundary conditions False is needed e.g. in plot routines to avoid showing incorrect bonds

  • tolerance (int) – tolerance (round decimal points) used for computing neighbor shells

  • id_list

  • cutoff (float/ None) – Upper bound of the distance to which the search must be done

  • cutoff_radius (float/ None) – Upper bound of the distance to which the search must be done

Returns

Neighbors instances with the neighbor indices, distances and vectors

Return type

pyiron.atomistics.structure.atoms.Neighbors

get_neighbors(num_neighbors=12, t_vec=True, include_boundary=True, exclude_self=True, tolerance=2, id_list=None, cutoff_radius=None, cutoff=None)[source]
Parameters
  • num_neighbors (int) – number of neighbors

  • t_vec (bool) – True: compute distance vectors (pbc are automatically taken into account)

  • include_boundary (bool) – True: search for neighbors assuming periodic boundary conditions False is needed e.g. in plot routines to avoid showing incorrect bonds

  • exclude_self (bool) – include central __atom (i.e. distance = 0)

  • tolerance (int) – tolerance (round decimal points) used for computing neighbor shells

  • id_list

  • cutoff (float/None) – Upper bound of the distance to which the search must be done - by default search for upto 100 neighbors unless num_neighbors is defined explicitly.

  • cutoff_radius (float/None) – Upper bound of the distance to which the search must be done - by default search for upto 100 neighbors unless num_neighbors is defined explicitly.

Returns

Neighbors instances with the neighbor indices, distances and vectors

Return type

pyiron.atomistics.structure.atoms.Neighbors

get_number_of_atoms()[source]

Returns:

get_number_of_degrees_of_freedom()[source]

Returns:

get_number_of_species()[source]

Returns:

get_number_species_atoms()[source]

Returns a dictionary with the species in the structure and the corresponding count in the structure

Returns

An ordered dictionary with the species and the corresponding count

Return type

collections.OrderedDict

get_parent_basis()[source]

Returns the basis with all user defined/special elements as the it’s parent

Returns

Structure without any user defined elements

Return type

pyiron.atomistics.structure.atoms.Atoms

get_parent_symbols()[source]

Returns the chemical symbols for all the atoms in the structure even for user defined elements

Returns

A list of chemical symbols

Return type

numpy.ndarray

get_pbc()[source]

Returns a boolean array of the periodic boundary conditions along the x, y and z axis respectively

Returns

Boolean array of length 3

Return type

numpy.ndarray

get_positions()[source]

Get positions. This function is for compatability with ASE

Returns

Positions in absolute coordinates

Return type

numpy.ndarray

get_primitive_cell(symprec=1e-05, angle_tolerance=- 1.0)[source]
Parameters
  • symprec

  • angle_tolerance

Returns:

get_scaled_positions(wrap=True)[source]

Returns the scaled/relative positions

Returns

The relative positions of the atoms in the supercell

Return type

numpy.ndarray

get_shell_matrix(shell=None, id_list=None, restraint_matrix=None, num_neighbors=100, tolerance=2)[source]
Parameters
  • neigh_list – user defined get_neighbors (recommended if atoms are displaced from the ideal positions)

  • id_list – cf. get_neighbors

  • radius – cf. get_neighbors

  • num_neighbors – cf. get_neighbors

  • tolerance – cf. get_neighbors

  • restraint_matrix – NxN matrix with True or False, where False will remove the entries. If an integer is given the sum of the chemical indices corresponding to the number will be set to True and the rest to False

Returns

NxN matrix with 1 for the pairs of atoms in the given shell

get_shell_radius(shell=1, id_list=None)[source]
Parameters
  • shell

  • id_list

Returns:

get_shells(id_list=None, max_shell=2, max_num_neighbors=100)[source]
Parameters
  • id_list

  • max_shell

  • max_num_neighbors

Returns:

get_spacegroup(symprec=1e-05, angle_tolerance=- 1.0)[source]
Parameters
  • symprec

  • angle_tolerance

Returns:

https://atztogo.github.io/spglib/python-spglib.html

get_species_objects()[source]

Returns:

get_species_symbols()[source]

Returns the symbols of the present species

Returns

List of the symbols of the species

Return type

numpy.ndarray

get_spherical_coordinates(x=None)[source]
Parameters

x (list/ndarray) – coordinates to transform. If not set, the positions in structure will be returned.

Returns

array in spherical coordinates

get_symmetry(use_magmoms=False, use_elements=True, symprec=1e-05, angle_tolerance=- 1.0)[source]
Parameters
  • use_magmoms

  • use_elements – True or False. If False, chemical elements will be ignored

  • symprec

  • angle_tolerance

Returns:

get_symmetry_dataset(symprec=1e-05, angle_tolerance=- 1.0)[source]
Parameters
  • symprec

  • angle_tolerance

Returns:

https://atztogo.github.io/spglib/python-spglib.html

get_tags()[source]

Returns the keys of the stored tags of the structure

Returns

Keys of the stored tags

Return type

dict_keys

get_volume(per_atom=False)[source]
Parameters

per_atom (bool) – True if volume per atom is to be returned

Returns

Volume in A**3

Return type

volume (float)

get_voronoi_volume()[source]

Returns:

group_points_by_symmetry(points)[source]

This function classifies the points into groups according to the box symmetry given by spglib.

Parameters

points – (np.array/list) nx3 array which contains positions

Returns: list of arrays containing geometrically equivalent positions

It is possible that the original points are not found in the returned list, as the positions outsie the box will be projected back to the box.

property info

This dictionary is merely used to be compatible with the ASE Atoms class.

Type

dict

new_array(name, a, dtype=None, shape=None)[source]

Adding a new array to the instance. This function is for the purpose of compatibility with the ASE package

Parameters
  • name (str) – Name of the array

  • a (list/numpy.ndarray) – The array to be added

  • dtype (type) – Data type of the array

  • shape (list/turple) – Shape of the array

numbers_to_elements(numbers)[source]

Convert atomic numbers in element objects (needed for compatibility with ASE)

Parameters

numbers (list) – List of Element Numbers (as Integers; default in ASE)

Returns

A list of elements as needed for pyiron

Return type

list

occupy_lattice(**qwargs)[source]

Replaces specified indices with a given species

property pbc

A list of boolean values which gives the periodic boundary consitions along the three axes. The default value is [True, True, True]

Type

list

plot3d(show_cell=True, show_axes=True, camera='orthographic', spacefill=True, particle_size=1.0, select_atoms=None, background='white', color_scheme=None, colors=None, scalar_field=None, scalar_start=None, scalar_end=None, scalar_cmap=None, vector_field=None, vector_color=None, magnetic_moments=False, custom_array=None, custom_3darray=None)[source]

Plot3d relies on NGLView to visualize atomic structures. Here, we construct a string in the “protein database” (“pdb”) format, then turn it into an NGLView “structure”. PDB is a white-space sensitive format, so the string snippets are carefully formatted.

The final widget is returned. If it is assigned to a variable, the visualization is suppressed until that variable is evaluated, and in the meantime more NGL operations can be applied to it to modify the visualization.

Parameters
  • show_cell (bool) – Whether or not to show the frame. (Default is True.)

  • show_axes (bool) – Whether or not to show xyz axes. (Default is True.)

  • camera (str) – ‘perspective’ or ‘orthographic’. (Default is ‘perspective’.)

  • spacefill (bool) – Whether to use a space-filling or ball-and-stick representation. (Default is True, use space-filling atoms.)

  • particle_size (float) – Size of the particles. (Default is 1.)

  • select_atoms (numpy.ndarray) – Indices of atoms to show, either as integers or a boolean array mask. (Default is None, show all atoms.)

  • background (str) – Background color. (Default is ‘white’.)

  • color_scheme (str) – NGLView color scheme to use. (Default is None, color by element.)

  • colors (numpy.ndarray) – A per-atom array of HTML color names or hex color codes to use for atomic colors. (Default is None, use coloring scheme.)

  • scalar_field (numpy.ndarray) – Color each atom according to the array value (Default is None, use coloring scheme.)

  • scalar_start (float) – The scalar value to be mapped onto the low end of the color map (lower values are clipped). (Default is None, use the minimum value in scalar_field.)

  • scalar_end (float) – The scalar value to be mapped onto the high end of the color map (higher values are clipped). (Default is None, use the maximum value in scalar_field.)

  • scalar_cmap (matplotlib.cm) – The colormap to use. (Default is None, giving a blue-red divergent map.)

  • vector_field (numpy.ndarray) – Add vectors (3 values) originating at each atom. (Default is None, no vectors.)

  • vector_color (numpy.ndarray) – Colors for the vectors (only available with vector_field). (Default is None, vectors are colored by their direction.)

  • magnetic_moments (bool) – Plot magnetic moments as ‘scalar_field’ or ‘vector_field’.

  • NGLView color schemes (Possible) – ” “, “picking”, “random”, “uniform”, “atomindex”, “residueindex”, “chainindex”, “modelindex”, “sstruc”, “element”, “resname”, “bfactor”, “hydrophobicity”, “value”, “volume”, “occupancy”

Returns

The NGLView widget itself, which can be operated on further or viewed as-is.

Return type

(nglview.NGLWidget)

Warning

  • Many features only work with space-filling atoms (e.g. coloring by a scalar field).

  • The colour interpretation of some hex codes is weird, e.g. ‘green’.

plot3d_ase(spacefill=True, show_cell=True, camera='perspective', particle_size=0.5, background='white', color_scheme='element', show_axes=True)[source]
Possible color schemes:

” “, “picking”, “random”, “uniform”, “atomindex”, “residueindex”, “chainindex”, “modelindex”, “sstruc”, “element”, “resname”, “bfactor”, “hydrophobicity”, “value”, “volume”, “occupancy”

Returns:

pop(i=- 1)[source]

Remove and return atom at index i (default last).

Parameters

i

Returns:

pos_xyz()[source]

Returns:

refine_cell(symprec=1e-05, angle_tolerance=- 1.0)[source]
Parameters
  • symprec

  • angle_tolerance

Returns:

https://atztogo.github.io/spglib/python-spglib.html

repeat(rep)[source]

Create new repeated atoms object.

The rep argument should be a sequence of three positive integers like (2,3,1) or a single integer (r) equivalent to (r,r,r).

repeat_points(points, rep, centered=False)[source]

Return points with repetition given according to periodic boundary conditions

Parameters
  • points (np.ndarray/list) – xyz vector or list/array of xyz vectors

  • rep (int/list/np.ndarray) – Repetition in each direction. If int is given, the same value is used for every direction

  • centered (bool) – Whether the original points should be in the center of repeated points.

Returns

(np.ndarray) repeated points

reset_absolute(is_absolute)[source]
rotate(vector, angle=None, center=0, 0, 0, rotate_cell=False, index_list=None)[source]

Rotate atoms based on a vector and an angle, or two vectors. This function is completely adopted from ASE code (https://wiki.fysik.dtu.dk/ase/_modules/ase/atoms.html#Atoms.rotate)

Parameters
  • rotate_cell

  • center

  • vector (list/numpy.ndarray/string) – Vector to rotate the atoms around. Vectors can be given as strings: ‘x’, ‘-x’, ‘y’, … .

  • angle (float/list) – Angle that the atoms is rotated around the vecor ‘v’. If an angle is not specified, the length of ‘v’ is used as the angle (default). The angle can also be a vector and then ‘v’ is rotated into ‘a’.

  • = [0 (center) – The center is kept fixed under the rotation. Use ‘COM’ to fix the center of mass, ‘COP’ to fix the center of positions or ‘COU’ to fix the center of cell.

  • 0 – The center is kept fixed under the rotation. Use ‘COM’ to fix the center of mass, ‘COP’ to fix the center of positions or ‘COU’ to fix the center of cell.

  • 0] – The center is kept fixed under the rotation. Use ‘COM’ to fix the center of mass, ‘COP’ to fix the center of positions or ‘COU’ to fix the center of cell.

  • = False (rotate_cell) – If true the cell is also rotated.

  • index_list (list/numpy.ndarray) – Indices of atoms to be rotated

Examples:

Rotate 90 degrees around the z-axis, so that the x-axis is rotated into the y-axis:

>>> atoms = Atoms('H', [[-0.1, 1.01, -0.5]], cell=[[1, 0, 0], [0, 1, 0], [0, 0, 4]], pbc=[1, 1, 0])
>>> a = (22./ 7.) / 2. # pi/2
>>> atoms.rotate('z', a)
>>> atoms.rotate((0, 0, 1), a)
>>> atoms.rotate('-z', -a)
>>> atoms.rotate((0, 0, a))
>>> atoms.rotate('x', 'y')
rotate_euler(center=0, 0, 0, phi=0.0, theta=0.0, psi=0.0)[source]

Rotate atoms via Euler angles.

See e.g http://mathworld.wolfram.com/EulerAngles.html for explanation.

Parameters:

center :

The point to rotate about. a sequence of length 3 with the coordinates, or ‘COM’ to select the center of mass, ‘COP’ to select center of positions or ‘COU’ to select center of cell.

phi :

The 1st rotation angle around the z axis (in radian)

theta :

Rotation around the x axis (in radian)

psi :

2nd rotation around the z axis (in radian)

scaled_pos_xyz()[source]

Returns:

property scaled_positions
select_index(el)[source]

Returns the indices of a given element in the structure

Parameters

el (str/atomistics.structures.periodic_table.ChemicalElement/list) – Element for which the indices should be returned

Returns

An array of indices of the atoms of the given element

Return type

numpy.ndarray

select_parent_index(el)[source]

Returns the indices of a given element in the structure ignoring user defined elements

Parameters

el (str/atomistics.structures.periodic_table.ChemicalElement) – Element for which the indices should be returned

Returns

An array of indices of the atoms of the given element

Return type

numpy.ndarray

set_absolute()[source]
set_array(name, a, dtype=None, shape=None)[source]

Update array. This function is for the purpose of compatibility with the ASE package

Parameters
  • name (str) – Name of the array

  • a (list/numpy.ndarray) – The array to be added

  • dtype (type) – Data type of the array

  • shape (list/turple) – Shape of the array

set_calculator(calc=None)[source]

Attach calculator object.

set_cell(cell, scale_atoms=False)[source]

Set unit cell vectors.

Parameters:

cell: 3x3 matrix or length 3 or 6 vector

Unit cell. A 3x3 matrix (the three unit cell vectors) or just three numbers for an orthorhombic cell. Another option is 6 numbers, which describes unit cell with lengths of unit cell vectors and with angles between them (in degrees), in following order: [len(a), len(b), len(c), angle(b,c), angle(a,c), angle(a,b)]. First vector will lie in x-direction, second in xy-plane, and the third one in z-positive subspace.

scale_atoms: bool

Fix atomic positions or move atoms with the unit cell? Default behavior is to not move the atoms (scale_atoms=False).

Examples:

Two equivalent ways to define an orthorhombic cell:

>>> atoms = Atoms('He')
>>> a, b, c = 7, 7.5, 8
>>> atoms.set_cell([a, b, c])
>>> atoms.set_cell([(a, 0, 0), (0, b, 0), (0, 0, c)])

FCC unit cell:

>>> atoms.set_cell([(0, b, b), (b, 0, b), (b, b, 0)])

Hexagonal unit cell:

>>> atoms.set_cell([a, a, c, 90, 90, 120])

Rhombohedral unit cell:

>>> alpha = 77
>>> atoms.set_cell([a, a, a, alpha, alpha, alpha])
set_constraint(constrain)[source]
set_initial_magnetic_moments(magmoms)[source]

Set array of initial magnetic moments.

Parameters

magmoms (numpy.array()) –

set_pbc(value)[source]

Sets the perioic boundary conditions on all three axis

Parameters

value (numpy.ndarray/list) – An array of bool type with length 3

set_positions(positions)[source]

Set positions. This function is for compatability with ASE

Parameters

positions (numpy.ndarray/list) – Positions in absolute coordinates

set_relative()[source]
set_repeat(vec)[source]
set_scaled_positions(scaled)[source]

Set positions relative to unit cell.

Parameters

scaled (numpy.ndarray/list) – The relative coordinates

set_species(value)[source]

Setting the species list

Parameters

value (list) – A list atomistics.structure.periodic_table.ChemicalElement instances

property species

A list of atomistics.structure.periodic_table.ChemicalElement instances

Type

list

symmetrize_vectors(vectors, force_update=False, use_magmoms=False, use_elements=True, symprec=1e-05, angle_tolerance=- 1.0)[source]

Symmetrization of natom x 3 vectors according to box symmetries

Parameters
  • vectors (ndarray/list) – natom x 3 array to symmetrize

  • force_update (bool) – whether to update the symmetry info

  • use_magmoms (bool) – cf. get_symmetry

  • use_elements (bool) – cf. get_symmetry

  • symprec (float) – cf. get_symmetry

  • angle_tolerance (float) – cf. get_symmetry

Returns

(np.ndarray) symmetrized vectors

to_hdf(hdf, group_name='structure')[source]

Save the object in a HDF5 file

Parameters
  • hdf (pyiron.base.generic.hdfio.FileHDFio) – HDF path to which the object is to be saved

  • group_name (str) – Group name with which the object should be stored. This same name should be used to retrieve the object

translate(displacement)[source]

Translate atomic positions.

The displacement argument can be a float, an xyz vector, or an nx3 array (where n is the number of atoms).

Parameters

displacement

Returns:

wrap(center=0.5, 0.5, 0.5, pbc=None, eps=1e-07)[source]

Wrap positions to unit cell.

Parameters:

center: three float

The positons in fractional coordinates that the new positions will be nearest possible to.

pbc: one or 3 bool

For each axis in the unit cell decides whether the positions will be moved along this axis. By default, the boundary conditions of the Atoms object will be used.

eps: float

Small number to prevent slightly negative coordinates from beeing wrapped.

See also the ase.utils.geometry.wrap_positions() function. Example:

>>> a = Atoms('H',
...           [[-0.1, 1.01, -0.5]],
...           cell=[[1, 0, 0], [0, 1, 0], [0, 0, 4]],
...           pbc=[1, 1, 0])
>>> a.wrap()
>>> a.positions
array([[ 0.9 ,  0.01, -0.5 ]])
write(filename, format=None, **kwargs)[source]

Write atoms object to a file.

see ase.io.write for formats. kwargs are passed to ase.io.write.

Parameters
  • filename

  • format

  • **kwargs

Returns:

class pyiron.atomistics.structure.atoms.CrystalStructure(*args, **kwargs)[source]

Bases: object

class pyiron.atomistics.structure.atoms.Neighbors[source]

Bases: object

Class for storage of the neighbor information for a given atom based on the KDtree algorithm

property distances
property indices
property shells
property vecs
pyiron.atomistics.structure.atoms.ase_to_pyiron(ase_obj)[source]
Parameters

ase_obj

Returns:

pyiron.atomistics.structure.atoms.default(data, dflt)[source]

Helper function for setting default values.

Parameters
  • data

  • dflt

Returns:

pyiron.atomistics.structure.atoms.ovito_to_pyiron(ovito_obj)[source]
Parameters

ovito_obj

Returns:

pyiron.atomistics.structure.atoms.pyiron_to_ase(pyiron_obj)[source]
pyiron.atomistics.structure.atoms.pyiron_to_ovito(atoms)[source]
Parameters

atoms

Returns:

pyiron.atomistics.structure.atoms.pyiron_to_pymatgen(pyiron_obj)[source]

Convert pyiron atoms object to pymatgen atoms object

Parameters

pyiron_obj – pyiron atoms object

Returns

pymatgen atoms object

pyiron.atomistics.structure.atoms.pymatgen_to_pyiron(pymatgen_obj)[source]

Convert pymatgen atoms object to pyiron atoms object (pymatgen->ASE->pyiron)

Parameters

pymatgen_obj – pymatgen atoms object

Returns

pyiron atoms object

pyiron.atomistics.structure.atoms.string2symbols(s)[source]

Convert string to list of chemical symbols.

Parameters

s

Returns:

pyiron.atomistics.structure.atoms.string2vector(v)[source]
Parameters

v

Returns:

pyiron.atomistics.structure.atoms.symbols2numbers(symbols)[source]
Parameters

symbols (list, str) –

Returns: