Geometry

class sisl.Geometry(xyz, atoms=None, sc=None, names=None)

Bases: sisl.supercell.SuperCellChild

Holds atomic information, coordinates, species, lattice vectors

The Geometry class holds information regarding atomic coordinates, the atomic species, the corresponding lattice-vectors.

It enables the interaction and conversion of atomic structures via simple routine methods.

All lengths are assumed to be in units of Angstrom, however, as long as units are kept same the exact units are irrespective.

>>> square = Geometry([[0.5, 0.5, 0.5]], Atom(1),
...                   sc=SuperCell([1, 1, 10], nsc=[3, 3, 1]))
>>> print(square)
Geometry{na: 1, no: 1,
 Atoms{species: 1,
  Atom{H, Z: 1, mass(au): 1.00794, maxR: -1.00000,
   Orbital{R: -1.00000, q0: 0.0}
  }: 1,
 },
 maxR: -1.00000,
 SuperCell{volume: 1.0000e+01, nsc: [3 3 1]}
}
na
xyz

atomic coordinates

Type

numpy.ndarray

atoms
orbitals
sc

the supercell describing the periodicity of the geometry

Type

SuperCell

no
n_s

total number of supercells in the supercell

Type

int

no_s

total number of orbitals in the geometry times number of supercells

Type

int

Parameters
  • xyz (array_like) – atomic coordinates xyz[i, :] is the atomic coordinate of the i’th atom.

  • atoms (array_like or Atoms) – atomic species retrieved from the PeriodicTable

  • sc (SuperCell) – the unit-cell describing the atoms in a periodic super-cell

Examples

An atomic cubic lattice of Hydrogen atoms

>>> xyz = [[0, 0, 0],
...        [1, 1, 1]]
>>> sc = SuperCell([2,2,2])
>>> g = Geometry(xyz, Atom('H'), sc)

The following estimates the lattice vectors from the atomic coordinates, although possible, it is not recommended to be used.

>>> xyz = [[0, 0, 0],
...        [1, 1, 1]]
>>> g = Geometry(xyz, Atom('H'))

See also

Atoms

contained atoms self.atoms

Atom

contained atoms are each an object of this

Attributes

__dict__

__doc__

__hash__

__module__

__weakref__

list of weak references to the object (if defined)

atom

deprecated atom

atoms

Atoms for the geometry (Atoms object)

cell

Returns the inherent SuperCell objects cell

firsto

The first orbital on the corresponding atom

fxyz

Returns geometry coordinates in fractional coordinates

icell

Returns the inherent SuperCell objects icell

isc_off

Returns the inherent SuperCell objects isc_off

lasto

The last orbital on the corresponding atom

mass

The mass of all atoms as an array

n_s

Returns the inherent SuperCell objects n_s

na

Number of atoms in geometry

na_s

Number of supercell atoms

names

The named index specifier

no

Number of orbitals

no_s

Number of supercell orbitals

nsc

Returns the inherent SuperCell objects nsc

orbitals

List of orbitals per atom

origo

Returns the inherent SuperCell objects origo

q0

Total initial charge in this geometry (sum of q0 in all atoms)

rcell

Returns the inherent SuperCell objects rcell

sc_off

Returns the inherent SuperCell objects sc_off

volume

Returns the inherent SuperCell objects vol

Methods

Rij(ia, ja)

Vector between atom ia and ja, atoms can be in super-cell indices

_(atoms)

_ArgumentParser_args_single()

Returns the options for Geometry.ArgumentParser in case they are the only options

_Geometry__currently_not_used_close_rec(xyz_ia)

Indices of atoms in a given supercell within a given radius from a given coordinate

_Geometry__init_sc(sc)

Initializes the supercell by calculating the size if not supplied

__add__(b)

Merge two geometries (or geometry and supercell)

__delattr__

Implement delattr(self, name).

__dir__

Default dir() implementation.

__eq__(other)

Return self==value.

__format__

Default object formatter.

__ge__

Return self>=value.

__getattribute__

Return getattr(self, name).

__getitem__(atoms)

Geometry coordinates (allows supercell indices)

__getstate__()

Returns the state of this object

__gt__

Return self>value.

__init__(xyz[, atoms, sc, names])

Initialize self.

__init_subclass__

This method is called when a class is subclassed.

__iter__()

An iterator over all atomic indices

__le__

Return self<=value.

__len__()

Number of atoms in geometry

__lt__

Return self<value.

__mul__(m)

Implement easy repeat function

__ne__(other)

Return self!=value.

__new__

Create and return a new object.

__plot__([axis, supercell, axes, atom_indices])

Plot the geometry in a specified matplotlib.Axes object.

__radd__(b)

Merge two geometries (or geometry and supercell)

__reduce__

Helper for pickle.

__reduce_ex__

Helper for pickle.

__repr__()

A simple, short string representation.

__rmul__(m)

Implement easy repeat function

__setattr__

Implement setattr(self, name, value).

__setitem__(atoms, value)

Specify geometry coordinates

__setstate__(d)

Re-create the state of this object

__sizeof__

Size of object in memory, in bytes.

__str__()

str of the object

__subclasshook__

Abstract classes can override this to customize issubclass().

_fill(non_filled[, dtype])

_fill_sc(supercell_index)

_sanitize_atoms(atoms)

Converts an atoms to index under given inputs

_sanitize_orbs(orbital)

Converts an orbital to index under given inputs

a2isc(atoms)

Returns super-cell index for a specific/list atom

a2o(atoms[, all])

Returns an orbital index of the first orbital of said atom.

a2sc(atoms)

Returns the super-cell offset for a specific atom

a2transpose(atoms1[, atoms2])

Transposes connections from atoms1 to atoms2 such that supercell connections are transposed

add(other[, offset])

Merge two geometries (or a Geometry and SuperCell) by adding the two atoms together

add_vacuum(vacuum, axis)

Add vacuum along the axis lattice vector

angle(atoms[, dir, ref, rad])

The angle between atom atoms and the direction dir, with possibility of a reference coordinate ref

append(other, axis[, offset])

Appends two structures along axis

area(ax0, ax1)

Calculate the area spanned by the two axis ax0 and ax1

as_primary(na_primary[, ret_super])

Try and reduce the geometry to the primary unit-cell comprising na_primary atoms

asc2uc(atoms[, unique])

Returns atoms from supercell indices to unit-cell indices, possibly removing dublicates

attach(atom, other, other_atom[, dist, axis])

Attaches another Geometry at the atom index with respect to other_atom using different methods.

auc2sc(atoms[, unique])

Returns atom from unit-cell indices to supercell indices, possibly removing dublicates

axyz([atoms, isc])

Return the atomic coordinates in the supercell of a given atom.

bond_correct(ia, atoms[, method])

Corrects the bond between ia and the atoms.

center([atoms, what])

Returns the center of the geometry

close(xyz_ia[, R, atoms, atoms_xyz, …])

Indices of atoms in the entire supercell within a given radius from a given coordinate

close_sc(xyz_ia[, isc, R, atoms, atoms_xyz, …])

Indices of atoms in a given supercell within a given radius from a given coordinate

copy()

A copy of the object.

cut(seps, axis[, seg, rtol, atol])

A subset of atoms from the geometry by cutting the geometry into seps parts along the direction axis.

distance([atoms, R, tol, method])

Calculate the distances for all atoms in shells of radius tol within max_R

equal(other[, R, tol])

Whether two geometries are the same (optional not check of the orbital radius)

fromASE(aseg)

Returns geometry from an ASE object.

iR([na, iR, R])

Return an integer number of maximum radii (self.maxR()) which holds approximately na atoms

insert(atom, geometry)

Inserts other atoms right before index

is_orthogonal()

Return true if all cell vectors are linearly independent

iter()

An iterator over all atomic indices

iter_block([iR, R, atoms, method])

Iterator for performance critical loops

iter_block_rand([iR, R, atoms])

Perform the random block-iteration by randomly selecting the next center of block

iter_block_shape([shape, iR, atoms])

Perform the grid block-iteration by looping a grid

iter_orbitals([atoms, local])

Returns an iterator over all atoms and their associated orbitals

iter_species([atoms])

Iterator over all atoms (or a subset) and species as a tuple in this geometry

maxR([all])

Maximum orbital range of the atoms

mirror(method[, atoms])

Mirrors the atomic coordinates about the plane

move(v[, atoms, cell])

Translates the geometry by v

o2a(io[, unique])

Atomic index corresponding to the orbital indicies.

o2isc(orbitals)

Returns the super-cell index for a specific orbital.

o2sc(orbitals)

Returns the super-cell offset for a specific orbital.

o2transpose(orb1[, orb2])

Transposes connections from orb1 to orb2 such that supercell connections are transposed

oRij(io, jo)

Vector between orbital io and jo, orbitals can be in super-cell indices

optimize_nsc([axis, R])

Optimize the number of supercell connections based on self.maxR()

orij(io, jo)

Distance between orbital io and jo, orbitals can be in super-cell indices

osc2uc(orb[, unique])

Returns orbitals from supercell indices to unit-cell indices, possibly removing dublicates

ouc2sc(orbitals[, unique])

Returns orbitals from unit-cell indices to supercell indices, possibly removing dublicates

overlap(other[, eps, offset, offset_other])

Calculate the overlapping indices between two geometries

prepend(other, axis[, offset])

Prepend two structures along axis

read(sile, *args, **kwargs)

Reads geometry from the Sile using Sile.read_geometry

reduce()

Remove all atoms not currently used in the self.atoms object

remove(atoms)

Remove atoms from the geometry.

reorder()

Reorders atoms according to first occurence in the geometry

repeat(reps, axis)

Create a repeated geometry

reverse([atoms])

Returns a reversed geometry

rij(ia, ja)

Distance between atom ia and ja, atoms can be in super-cell indices

rotate(angle, v[, origo, atoms, only, rad])

Rotate geometry around vector and return a new geometry

rotate_miller(m, v)

Align Miller direction along v

sc2uc(atoms[, unique])

Returns atoms from supercell indices to unit-cell indices, possibly removing dublicates

sc_index(*args, **kwargs)

Call local SuperCell object sc_index function

scale(scale)

Scale coordinates and unit-cell to get a new geometry with proper scaling

set_nsc(*args, **kwargs)

Set the number of super-cells in the SuperCell object

set_sc(sc)

Overwrites the local supercell

set_supercell(sc)

Overwrites the local supercell

sort(**kwargs)

Sort atoms in a nested fashion according to various criteria

sparserij([dtype, na_iR, method])

Return the sparse matrix with all distances in the matrix

sub(atoms[, cell])

Create a new Geometry with a subset of this Geometry

swap(atoms_a, atoms_b)

Swap a set of atoms in the geometry and return a new one

swapaxes(axis_a, axis_b[, swap])

Swap the axis for the atomic coordinates and the cell vectors

tile(reps, axis)

Tile the geometry to create a bigger one

toASE()

Returns the geometry as an ASE Atoms object

translate(v[, atoms, cell])

Translates the geometry by v

uc2sc(atoms[, unique])

Returns atom from unit-cell indices to supercell indices, possibly removing dublicates

within(shapes[, atoms, atoms_xyz, ret_xyz, …])

Indices of atoms in the entire supercell within a given shape from a given coordinate

within_inf(sc[, periodic, tol, origo])

Find all atoms within a provided supercell

within_sc(shapes[, isc, atoms, atoms_xyz, …])

Indices of atoms in a given supercell within a given shape from a given coordinate

write(sile, *args, **kwargs)

Writes geometry to the Sile using sile.write_geometry

Rij(ia, ja)[source]

Vector between atom ia and ja, atoms can be in super-cell indices

Returns the vector between two atoms:

\[R_{ij} = r_j - r_i\]
Parameters
  • ia (int or array_like) – atomic index of first atom

  • ja (int or array_like) – atomic indices

a2isc(atoms)[source]

Returns super-cell index for a specific/list atom

Returns a vector of 3 numbers with integers.

a2o(atoms, all=False)[source]

Returns an orbital index of the first orbital of said atom. This is particularly handy if you want to create TB models with more than one orbital per atom.

Note that this will preserve the super-cell offsets.

Parameters
  • atoms (array_like) – Atomic indices

  • all (bool, optional) – False, return only the first orbital corresponding to the atom, True, returns list of the full atom

a2sc(atoms)[source]

Returns the super-cell offset for a specific atom

a2transpose(atoms1, atoms2=None)[source]

Transposes connections from atoms1 to atoms2 such that supercell connections are transposed

When handling supercell indices it is useful to get the transposed connection. I.e. if you have a connection from site i (in unit cell indices) to site j (in supercell indices) it may be useful to get the equivalent supercell connection such for site j (in unit cell indices) to site i (in supercell indices) such that they correspond to the transposed coupling.

Note that since this transposes couplings the indices returned are always expanded to the full length if either of the inputs are a single index.

Examples

>>> gr = geom.graphene()
>>> atoms = gr.close(0, 1.5)
>>> atoms
array([0, 1, 5, 9], dtype=int32)
>>> gr.a2transpose(0, atoms)
(array([0, 1, 1, 1], dtype=int32), array([ 0,  0, 14, 10], dtype=int32))
Parameters
  • atoms1 (array_like) – atomic indices must have same length as atoms2 or length 1

  • atoms2 (array_like, optional) – atomic indices must have same length as atoms1 or length 1. If not present then only atoms1 will be returned in transposed indices.

Returns

  • atoms2 (array_like) – transposed indices for atoms2 (only returned if atoms2 is not None)

  • atoms1 (array_like) – transposed indices for atoms1

add(other, offset=0, 0, 0)[source]

Merge two geometries (or a Geometry and SuperCell) by adding the two atoms together

If other is a Geometry only the atoms gets added, to also add the supercell vectors simply do geom.add(other).add(other.sc).

Parameters
  • other (Geometry or SuperCell) – Other geometry class which is added

  • offset ((3,), optional) – offset in geometry of other when adding the atoms. Only if other is of instance Geometry.

See also

append

appending geometries

prepend

prending geometries

attach

attach a geometry

insert

insert a geometry

add_vacuum(vacuum, axis)

Add vacuum along the axis lattice vector

Parameters
  • vacuum (float) – amount of vacuum added, in Ang

  • axis (int) – the lattice vector to add vacuum along

angle(atoms, dir=1.0, 0, 0, ref=None, rad=False)[source]

The angle between atom atoms and the direction dir, with possibility of a reference coordinate ref

The calculated angle can be written as this

\[\alpha = \arccos \frac{(\mathrm{atom} - \mathrm{ref})\cdot \mathrm{dir}} {|\mathrm{atom}-\mathrm{ref}||\mathrm{dir}|}\]

and thus lies in the interval \([0 ; \pi]\) as one cannot distinguish orientation without additional vectors.

Parameters
  • atoms (int or array_like) – indices/boolean of all atoms where angles should be calculated on

  • dir (str, int or array_like, optional) – the direction from which the angle is calculated from, default to x

  • ref (int or array_like, optional) – the reference point from which the vectors are drawn, default to origo

  • rad (bool, optional) – whether the returned value is in radians

append(other, axis, offset='none')[source]

Appends two structures along axis

This will automatically add the self.cell[axis,:] to all atomic coordiates in the other structure before appending.

The basic algorithm is this:

>>> oxa = other.xyz + self.cell[axis,:][None,:]
>>> self.xyz = np.append(self.xyz,oxa)
>>> self.cell[axis,:] += other.cell[axis,:]

NOTE: The cell appended is only in the axis that is appended, which means that the other cell directions need not conform.

Parameters
  • other (Geometry or SuperCell) – Other geometry class which needs to be appended If a SuperCell only the super cell will be extended

  • axis (int) – Cell direction to which the other geometry should be appended.

  • offset ({'none', 'min', (3,)}) – By default appending two structures will simply use the coordinates, as is. With ‘min’, the routine will shift both the structures along the cell axis of self such that they coincide at the first atom, lastly one may use a specified offset to manually select how other is displaced. NOTE: That self.cell[axis, :] will be added to offset always.

See also

add

add geometries

prepend

prending geometries

attach

attach a geometry

insert

insert a geometry

area(ax0, ax1)

Calculate the area spanned by the two axis ax0 and ax1

as_primary(na_primary, ret_super=False)[source]

Try and reduce the geometry to the primary unit-cell comprising na_primary atoms

This will basically try and find the tiling/repetitions required for the geometry to only have na_primary atoms in the unit cell.

Parameters
  • na_primary (int) – number of atoms in the primary unit cell

  • ret_super (bool, optional) – also return the number of supercells used in each direction

Returns

  • Geometry – the primary unit cell

  • SuperCell – the tiled supercell numbers used to find the primary unit cell (only if ret_super is true)

Raises

SislError – If the algorithm fails.

asc2uc(atoms, unique=False)

Returns atoms from supercell indices to unit-cell indices, possibly removing dublicates

Parameters
  • atoms (array_like or int) – the atomic supercell indices to be converted to unit-cell indices

  • unique (bool, optional) – If True the returned indices are unique and sorted.

property atom

deprecated atom

property atoms

Atoms for the geometry (Atoms object)

attach(atom, other, other_atom, dist='calc', axis=None)[source]

Attaches another Geometry at the atom index with respect to other_atom using different methods.

The attached geometry will be inserted at the end of the geometry via add.

Parameters
  • atom (int) – atomic index which is the base position of the attachment. The distance between atom and other_atom is dist.

  • other (Geometry) – the other Geometry to attach at the given point. In this case dist from atom.

  • other_atom (int) – the index of the atom in other that is inserted at atom.

  • dist (array_like or float or str, optional) – the distance (in Ang) between the attached coordinates. If dist is array_like it should be the vector between the atoms; if dist is float the argument axis is required and the vector will be calculated along the corresponding latticevector; else if dist is str this will correspond to the method argument of the Atom.radius class of the two atoms. Here axis is also required.

  • axis (int) – specify the direction of the lattice vectors used. Not used if dist is an array-like argument.

auc2sc(atoms, unique=False)

Returns atom from unit-cell indices to supercell indices, possibly removing dublicates

Parameters
  • atoms (array_like or int) – the atomic unit-cell indices to be converted to supercell indices

  • unique (bool, optional) – If True the returned indices are unique and sorted.

axyz(atoms=None, isc=None)[source]

Return the atomic coordinates in the supercell of a given atom.

The Geometry[...] slicing is calling this function with appropriate options.

Parameters
  • atoms (int or array_like) – atom(s) from which we should return the coordinates, the atomic indices may be in supercell format.

  • isc (array_like, optional) – Returns the atomic coordinates shifted according to the integer parts of the cell. Defaults to the unit-cell

Examples

>>> geom = Geometry([[0, 0, 0], [0.5, 0, 0]], sc=1.)
>>> print(geom.axyz(isc=[1,0,0])) 
[[1.   0.   0. ]
 [1.5  0.   0. ]]
>>> geom = Geometry([[0, 0, 0], [0.5, 0, 0]], sc=1.)
>>> print(geom.axyz(0)) 
[0.  0.  0.]
bond_correct(ia, atoms, method='calc')[source]

Corrects the bond between ia and the atoms.

Corrects the bond-length between atom ia and atoms in such a way that the atomic radius is preserved. I.e. the sum of the bond-lengths minimizes the distance matrix.

Only atom ia is moved.

Parameters
  • ia (int) – The atom to be displaced according to the atomic radius

  • atoms (array_like or int) – The atom(s) from which the radius should be reduced.

  • method (str, float, optional) – If str will use that as lookup in Atom.radius. Else it will be the new bond-length.

property cell

Returns the inherent SuperCell objects cell

center(atoms=None, what='xyz')[source]

Returns the center of the geometry

By specifying what one can control whether it should be:

  • xyz|position: Center of coordinates (default)

  • mm(xyz): Center of minimum/maximum of coordinates

  • mass: Center of mass

  • cell: Center of cell

Parameters
  • atoms (array_like) – list of atomic indices to find center of

  • what ({'xyz', 'mm(xyz)', 'mass', 'cell'}) – determine whether center should be of ‘cell’, mass-centered (‘mass’), center of minimum/maximum position of atoms or absolute center of the positions.

close(xyz_ia, R=None, atoms=None, atoms_xyz=None, ret_xyz=False, ret_rij=False)[source]

Indices of atoms in the entire supercell within a given radius from a given coordinate

This heavily relies on the close_sc method.

Note that if a connection is made in a neighbouring super-cell then the atomic index is shifted by the super-cell index times number of atoms. This allows one to decipher super-cell atoms from unit-cell atoms.

Parameters
  • xyz_ia (coordinate/index) – Either a point in space or an index of an atom. If an index is passed it is the equivalent of passing the atomic coordinate close_sc(self.xyz[xyz_ia,:]).

  • R ((None), float/tuple of float) –

    The radii parameter to where the atomic connections are found. If R is an array it will return the indices: in the ranges:

    >>> ( x <= R[0] , R[0] < x <= R[1], R[1] < x <= R[2] )
    

    If a single float it will return:

    >>> x <= R
    

  • atoms (array_like, optional) – List of indices for atoms that are to be considered

  • atoms_xyz (array_like, optional) – The atomic coordinates of the equivalent atoms variable (atoms must also be passed)

  • ret_xyz (bool, optional) – If true this method will return the coordinates for each of the couplings.

  • ret_rij (bool, optional) – If true this method will return the distances from the xyz_ia for each of the couplings.

Returns

  • index – indices of atoms (in supercell indices) within the shells of radius R

  • xyz – atomic coordinates of the indexed atoms (only for true ret_xyz)

  • rij – distance of the indexed atoms to the center coordinate (only for true ret_rij)

close_sc(xyz_ia, isc=0, 0, 0, R=None, atoms=None, atoms_xyz=None, ret_xyz=False, ret_rij=False)[source]

Indices of atoms in a given supercell within a given radius from a given coordinate

This returns a set of atomic indices which are within a sphere of radius R.

If R is a tuple/list/array it will return the indices: in the ranges:

>>> ( x <= R[0] , R[0] < x <= R[1], R[1] < x <= R[2] )
Parameters
  • xyz_ia (array_like of floats or int) – Either a point in space or an index of an atom. If an index is passed it is the equivalent of passing the atomic coordinate close_sc(self.xyz[xyz_ia,:]).

  • isc (array_like, optional) – The super-cell which the coordinates are checked in.

  • R (float or array_like, optional) – The radii parameter to where the atomic connections are found. If R is an array it will return the indices: in the ranges ( x <= R[0] , R[0] < x <= R[1], R[1] < x <= R[2] ). If a single float it will return x <= R.

  • atoms (array_like of int, optional) – List of atoms that will be considered. This can be used to only take out a certain atoms.

  • atoms_xyz (array_like of float, optional) – The atomic coordinates of the equivalent atoms variable (atoms must also be passed)

  • ret_xyz (bool, optional) – If True this method will return the coordinates for each of the couplings.

  • ret_rij (bool, optional) – If True this method will return the distance for each of the couplings.

Returns

  • index – indices of atoms (in supercell indices) within the shells of radius R

  • xyz – atomic coordinates of the indexed atoms (only for true ret_xyz)

  • rij – distance of the indexed atoms to the center coordinate (only for true ret_rij)

copy()[source]

A copy of the object.

cut(seps, axis, seg=0, rtol=0.0001, atol=0.0001)[source]

A subset of atoms from the geometry by cutting the geometry into seps parts along the direction axis.

This will effectively change the unit-cell in the axis as-well as removing self.na/seps atoms. It requires that self.na % seps == 0.

REMARK: You need to ensure that all atoms within the first cut out region are within the primary unit-cell.

Doing geom.cut(2, 1).tile(2, 1), could for symmetric setups, be equivalent to a no-op operation. A UserWarning will be issued if this is not the case.

This method may be regarded as the opposite of tile.

Parameters
  • seps (int) – number of times the structure will be cut.

  • axis (int) – the axis that will be cut

  • seg (int, optional) – returns the i’th segment of the cut structure Currently the atomic coordinates are not translated, this may change in the future.

  • rtol ((tolerance for checking tiling, see numpy.allclose)) –

  • atol ((tolerance for checking tiling, see numpy.allclose)) –

Examples

>>> g = sisl.geom.graphene()
>>> gxyz = g.tile(4, 0).tile(3, 1).tile(2, 2)
>>> G = gxyz.cut(2, 2).cut(3, 1).cut(4, 0)
>>> np.allclose(g.xyz, G.xyz)
True

See also

tile

opposite method of this

distance(atoms=None, R=None, tol=0.1, method='average')[source]

Calculate the distances for all atoms in shells of radius tol within max_R

Parameters
  • atoms (int or array_like, optional) – only create list of distances from the given atoms, default to all atoms

  • R (float, optional) – the maximum radius to consider, default to self.maxR(). To retrieve all distances for atoms within the supercell structure you can pass numpy.inf.

  • tol (float or array_like, optional) –

    the tolerance for grouping a set of atoms. This parameter sets the shell radius for each shell. I.e. the returned distances between two shells will be maximally 2*tol, but only if atoms are within two consecutive lists. If this is a list, the shells will be of unequal size.

    The first shell size will be tol * .5 or tol[0] * .5 if tol is a list.

  • method ({'average', 'mode', '<numpy.func>', func}) – How the distance in each shell is determined. A list of distances within each shell is gathered and the equivalent method will be used to extract a single quantity from the list of distances in the shell. If ‘mode’ is chosen it will use scipy.stats.mode. If a string is given it will correspond to getattr(numpy, method), while any callable function may be passed. The passed function will only be passed a list of unsorted distances that needs to be processed.

Examples

>>> geom = Geometry([0]*3, Atom(1, R=1.), sc=SuperCell(1., nsc=[5, 5, 1]))
>>> geom.distance() # use geom.maxR() 
array([1.])
>>> geom.distance(tol=[0.5, 0.4, 0.3, 0.2])
array([1.])
>>> geom.distance(R=2, tol=[0.5, 0.4, 0.3, 0.2]) 
array([1.        ,  1.41421356,  2.        ])
>>> geom.distance(R=2, tol=[0.5, 0.7]) # the R = 1 and R = 2 ** .5 gets averaged 
array([1.20710678,  2.        ])
Returns

an array of positive numbers yielding the distances from the atoms in reduced form

Return type

numpy.ndarray

See also

sparserij

return a sparse matrix will all distances between atoms

equal(other, R=True, tol=0.0001)[source]

Whether two geometries are the same (optional not check of the orbital radius)

Parameters
  • other (Geometry) – the other Geometry to check against

  • R (bool, optional) – if True also check if the orbital radii are the same (see Atom.equal)

  • tol (float, optional) – tolerance for checking the atomic coordinates

property firsto

The first orbital on the corresponding atom

classmethod fromASE(aseg)[source]

Returns geometry from an ASE object.

Parameters

aseg (ASE Atoms object which contains the following routines:) – get_atomic_numbers, get_positions, get_cell. From those methods a Geometry object will be created.

property fxyz

Returns geometry coordinates in fractional coordinates

iR(na=1000, iR=20, R=None)[source]

Return an integer number of maximum radii (self.maxR()) which holds approximately na atoms

Parameters
  • na (int, optional) – number of atoms within the radius

  • iR (int, optional) – initial iR value, which the sphere is estitametd from

  • R (float, optional) – the value used for atomic range (defaults to self.maxR())

Returns

number of radius needed to contain na atoms. Minimally 2 will be returned.

Return type

int

property icell

Returns the inherent SuperCell objects icell

insert(atom, geometry)[source]

Inserts other atoms right before index

We insert the geometry Geometry before atom. Note that this will not change the unit cell.

Parameters
  • atom (int) – the atomic index at which the other geometry is inserted

  • geometry (Geometry) – the other geometry to be inserted

See also

add

add geometries

append

appending geometries

prepend

prending geometries

attach

attach a geometry

is_orthogonal()

Return true if all cell vectors are linearly independent

property isc_off

Returns the inherent SuperCell objects isc_off

iter()[source]

An iterator over all atomic indices

This iterator is the same as:

>>> for ia in range(len(self)):
...    <do something>

or equivalently

>>> for ia in self:
...    <do something>

See also

iter_species

iterate across indices and atomic species

iter_orbitals

iterate across atomic indices and orbital indices

iter_block(iR=20, R=None, atoms=None, method='rand')[source]

Iterator for performance critical loops

NOTE: This requires that R has been set correctly as the maximum interaction range.

I.e. the loop would look like this:

>>> for ias, idxs in self.iter_block():
...    for ia in ias:
...        idx_a = self.close(ia, R = R, idx = idxs)

This iterator is intended for systems with more than 1000 atoms.

Remark that the iterator used is non-deterministic, i.e. any two iterators need not return the same atoms in any way.

Parameters
  • iR (int, optional) – the number of R ranges taken into account when doing the iterator

  • R (float, optional) – enables overwriting the local R quantity. Defaults to self.maxR()

  • atoms (array_like, optional) – enables only effectively looping a subset of the full geometry

  • method ({'rand', 'sphere', 'cube'}) –

    select the method by which the block iteration is performed. Possible values are:

    rand: a spherical object is constructed with a random center according to the internal atoms sphere: a spherical equispaced shape is constructed and looped cube: a cube shape is constructed and looped

Returns

  • numpy.ndarray – current list of atoms currently searched

  • numpy.ndarray – atoms that needs searching

iter_block_rand(iR=20, R=None, atoms=None)[source]

Perform the random block-iteration by randomly selecting the next center of block

iter_block_shape(shape=None, iR=20, atoms=None)[source]

Perform the grid block-iteration by looping a grid

iter_orbitals(atoms=None, local=True)[source]

Returns an iterator over all atoms and their associated orbitals

>>> for ia, io in self.iter_orbitals():

with ia being the atomic index, io the associated orbital index on atom ia. Note that io will start from 0.

Parameters
  • atoms (int or array_like, optional) – only loop on the given atoms, default to all atoms

  • local (bool, optional) – whether the orbital index is the global index, or the local index relative to the atom it resides on.

Yields
  • ia – atomic index

  • io – orbital index

See also

iter

iterate over atomic indices

iter_species

iterate across indices and atomic species

iter_species(atoms=None)[source]

Iterator over all atoms (or a subset) and species as a tuple in this geometry

>>> for ia, a, idx_specie in self.iter_species():
...     isinstance(ia, int) == True
...     isinstance(a, Atom) == True
...     isinstance(idx_specie, int) == True

with ia being the atomic index, a the Atom object, idx_specie is the index of the specie

Parameters

atoms (int or array_like, optional) – only loop on the given atoms, default to all atoms

See also

iter

iterate over atomic indices

iter_orbitals

iterate across atomic indices and orbital indices

property lasto

The last orbital on the corresponding atom

property mass

The mass of all atoms as an array

maxR(all=False)[source]

Maximum orbital range of the atoms

mirror(method, atoms=None)[source]

Mirrors the atomic coordinates about the plane

This will typically move the atomic coordinates outside of the unit-cell. This method should be used with care.

Parameters
  • method ({'xy'/'z', .., 'ab', .., v}) – mirror the structure about a Cartesian direction (x, y, z), plane (xy, xz, yz) or about user defined vectors (v). A vector may also be specified by 'ab' which is the vector normal to the plane spanned by the first and second lattice vector. or user defined vector (v)

  • atoms (array_like, optional) – only mirror a subset of atoms

move(v, atoms=None, cell=False)[source]

Translates the geometry by v

One can translate a subset of the atoms by supplying atoms.

Returns a copy of the structure translated by v.

Parameters
  • v (array_like) – the vector to displace all atomic coordinates

  • atoms (int or array_like, optional) – only displace the given atomic indices, if not specified, all atoms will be displaced

  • cell (bool, optional) – If True the supercell also gets enlarged by the vector

property n_s

Returns the inherent SuperCell objects n_s

property na

Number of atoms in geometry

property na_s

Number of supercell atoms

property names

The named index specifier

property no

Number of orbitals

property no_s

Number of supercell orbitals

property nsc

Returns the inherent SuperCell objects nsc

o2a(io, unique=False)[source]

Atomic index corresponding to the orbital indicies.

This is a particurlaly slow algorithm due to for-loops.

Note that this will preserve the super-cell offsets.

Parameters
  • io (array_like) – List of indices to return the atoms for

  • unique (bool, optional) – If True only return the unique atoms.

o2isc(orbitals)[source]

Returns the super-cell index for a specific orbital.

Returns a vector of 3 numbers with integers.

o2sc(orbitals)[source]

Returns the super-cell offset for a specific orbital.

o2transpose(orb1, orb2=None)[source]

Transposes connections from orb1 to orb2 such that supercell connections are transposed

When handling supercell indices it is useful to get the transposed connection. I.e. if you have a connection from site i (in unit cell indices) to site j (in supercell indices) it may be useful to get the equivalent supercell connection such for site j (in unit cell indices) to site i (in supercell indices) such that they correspond to the transposed coupling.

Note that since this transposes couplings the indices returned are always expanded to the full length if either of the inputs are a single index.

Examples

>>> gr = geom.graphene() # one orbital per site
>>> atoms = gr.close(0, 1.5)
>>> atoms
array([0, 1, 5, 9], dtype=int32)
>>> gr.o2transpose(0, atoms)
(array([0, 1, 1, 1], dtype=int32), array([ 0,  0, 14, 10], dtype=int32))
Parameters
  • orb1 (array_like) – orbital indices must have same length as orb2 or length 1

  • orb2 (array_like, optional) – orbital indices must have same length as orb1 or length 1. If not present then only orb1 will be returned in transposed indices.

Returns

  • orb2 (array_like) – transposed indices for orb2 (only returned if orb2 is not None)

  • orb1 (array_like) – transposed indices for orb1

oRij(io, jo)[source]

Vector between orbital io and jo, orbitals can be in super-cell indices

Returns the vector between two orbitals:

\[R_{ij} = r_j - r_i\]
Parameters
  • io (int or array_like) – orbital index of first orbital

  • jo (int or array_like) – orbital indices

optimize_nsc(axis=None, R=None)[source]

Optimize the number of supercell connections based on self.maxR()

After this routine the number of supercells may not necessarily be the same.

This is an in-place operation.

Parameters
  • axis (int or array_like, optional) – only optimize the specified axis (default to all)

  • R (float, optional) – the maximum connection radius for each atom

property orbitals

List of orbitals per atom

property origo

Returns the inherent SuperCell objects origo

orij(io, jo)[source]

Distance between orbital io and jo, orbitals can be in super-cell indices

Returns the distance between two orbitals:

\[r_{ij} = |r_j - r_i|\]
Parameters
  • io (int or array_like) – orbital index of first orbital

  • jo (int or array_like) – orbital indices

osc2uc(orb, unique=False)[source]

Returns orbitals from supercell indices to unit-cell indices, possibly removing dublicates

Parameters
  • orb (array_like or int) – the orbital supercell indices to be converted to unit-cell indices

  • unique (bool, optional) – If True the returned indices are unique and sorted.

ouc2sc(orbitals, unique=False)[source]

Returns orbitals from unit-cell indices to supercell indices, possibly removing dublicates

Parameters
  • orbitals (array_like or int) – the orbital unit-cell indices to be converted to supercell indices

  • unique (bool, optional) – If True the returned indices are unique and sorted.

overlap(other, eps=0.1, offset=0.0, 0.0, 0.0, offset_other=0.0, 0.0, 0.0)[source]

Calculate the overlapping indices between two geometries

Find equivalent atoms (in the primary unit-cell only) in two geometries. This routine finds which atoms have the same atomic positions in self and other.

Note that this will return duplicate overlapping atoms if one atoms lies within eps of more than 1 atom in other.

Parameters
  • other (Geometry) – Geometry to compare with self

  • eps (float, optional) – atoms within this distance will be considered equivalent

  • offset (list of float, optional) – offset for self.xyz before comparing

  • offset_other (list of float, optional) – offset for other.xyz before comparing

Examples

>>> gr22 = sisl.geom.graphene().tile(2, 0).tile(2, 1)
>>> gr44 = gr22.tile(2, 0).tile(2, 1)
>>> offset = np.array([0.2, 0.4, 0.4])
>>> gr22 = gr22.translate(offset)
>>> idx = np.arange(len(gr22))
>>> np.random.shuffle(idx)
>>> gr22 = gr22.sub(idx)
>>> idx22, idx44 = gr22.overlap(gr44, offset=-offset)
>>> assert idx22 == np.arange(len(gr22))
>>> assert idx44 == idx
Returns

  • idx_self (numpy.ndarray of int) – indices in self that are equivalent with idx_other

  • idx_other (numpy.ndarray of int) – indices in other that are equivalent with idx_self

prepend(other, axis, offset='none')[source]

Prepend two structures along axis

This will automatically add the self.cell[axis,:] to all atomic coordiates in the other structure before appending.

The basic algorithm is this:

>>> oxa = other.xyz
>>> self.xyz = np.append(oxa, self.xyz + other.cell[axis,:][None,:])
>>> self.cell[axis,:] += other.cell[axis,:]

NOTE: The cell prepended is only in the axis that is prependend, which means that the other cell directions need not conform.

Parameters
  • other (Geometry or SuperCell) – Other geometry class which needs to be prepended If a SuperCell only the super cell will be extended

  • axis (int) – Cell direction to which the other geometry should be prepended

  • offset ({'none', 'min', (3,)}) – By default appending two structures will simply use the coordinates, as is. With ‘min’, the routine will shift both the structures along the cell axis of other such that they coincide at the first atom, lastly one may use a specified offset to manually select how self is displaced. NOTE: That other.cell[axis, :] will be added to offset always.

See also

add

add geometries

append

appending geometries

attach

attach a geometry

insert

insert a geometry

property q0

Total initial charge in this geometry (sum of q0 in all atoms)

property rcell

Returns the inherent SuperCell objects rcell

static read(sile, *args, **kwargs)[source]

Reads geometry from the Sile using Sile.read_geometry

Parameters

sile (Sile, str or pathlib.Path) – a Sile object which will be used to read the geometry if it is a string it will create a new sile using get_sile.

See also

write

writes a Geometry to a given Sile/file

reduce()[source]

Remove all atoms not currently used in the self.atoms object

Notes

This is an in-place operation.

remove(atoms)[source]

Remove atoms from the geometry.

Indices passed MUST be unique.

Negative indices are wrapped and thus works.

Parameters

atoms (int or array_like) – indices/boolean of all atoms to be removed

See also

sub

the negative of this routine, i.e. retain a subset of atoms

reorder()[source]

Reorders atoms according to first occurence in the geometry

Notes

This is an in-place operation.

repeat(reps, axis)[source]

Create a repeated geometry

The atomic indices are NOT retained from the base structure.

The expansion of the atoms are basically performed using this algorithm:

>>> ja = 0
>>> for ia in range(self.na):
...     for id,r in args:
...        for i in range(r):
...           ja = ia + cell[id,:] * i

This method allows to utilise Bloch’s theorem when creating Hamiltonian parameter sets for TBtrans.

For geometries with a single atom this routine returns the same as tile.

Tiling and repeating a geometry will result in the same geometry. The only difference between the two is the final ordering of the atoms.

Parameters
  • reps (int) – number of repetitions

  • axis (int) – direction of repetition, 0, 1, 2 according to the cell-direction

Examples

>>> geom = Geometry([[0, 0, 0], [0.5, 0, 0]], sc=1)
>>> g = geom.repeat(2,axis=0)
>>> print(g.xyz) 
[[0.   0.   0. ]
 [1.   0.   0. ]
 [0.5  0.   0. ]
 [1.5  0.   0. ]]
>>> g = geom.repeat(2,0).repeat(2,1)
>>> print(g.xyz) 
[[0.   0.   0. ]
 [0.   1.   0. ]
 [1.   0.   0. ]
 [1.   1.   0. ]
 [0.5  0.   0. ]
 [0.5  1.   0. ]
 [1.5  0.   0. ]
 [1.5  1.   0. ]]

See also

tile

equivalent but different ordering of final structure

reverse(atoms=None)[source]

Returns a reversed geometry

Also enables reversing a subset of the atoms.

Parameters

atoms (int or array_like, optional) – only reverse the given atomic indices, if not specified, all atoms will be reversed

rij(ia, ja)[source]

Distance between atom ia and ja, atoms can be in super-cell indices

Returns the distance between two atoms:

\[r_{ij} = |r_j - r_i|\]
Parameters
  • ia (int or array_like) – atomic index of first atom

  • ja (int or array_like) – atomic indices

rotate(angle, v, origo=None, atoms=None, only='abc+xyz', rad=False)[source]

Rotate geometry around vector and return a new geometry

Per default will the entire geometry be rotated, such that everything is aligned as before rotation.

However, by supplying only = 'abc|xyz' one can designate which part of the geometry that will be rotated.

Parameters
  • angle (float) – the angle in degrees to rotate the geometry. Set the rad argument to use radians.

  • v (array_like) – the normal vector to the rotated plane, i.e. v = [1,0,0] will rotate the yz plane

  • origo (int or array_like, optional) – the origin of rotation. Anything but [0, 0, 0] is equivalent to a self.move(-origo).rotate(…).move(origo). If this is an int it corresponds to the atomic index.

  • atoms (int or array_like, optional) – only rotate the given atomic indices, if not specified, all atoms will be rotated.

  • only ({'abc+xyz', 'xyz', 'abc'}) – which coordinate subject should be rotated, if abc is in this string the cell will be rotated if xyz is in this string the coordinates will be rotated

  • rad (bool, optional) – if True the angle is provided in radians (rather than degrees)

See also

Quaternion

class to rotate

rotate_miller(m, v)[source]

Align Miller direction along v

Rotate geometry and cell such that the Miller direction points along the Cartesian vector v.

sc2uc(atoms, unique=False)[source]

Returns atoms from supercell indices to unit-cell indices, possibly removing dublicates

Parameters
  • atoms (array_like or int) – the atomic supercell indices to be converted to unit-cell indices

  • unique (bool, optional) – If True the returned indices are unique and sorted.

sc_index(*args, **kwargs)

Call local SuperCell object sc_index function

property sc_off

Returns the inherent SuperCell objects sc_off

scale(scale)[source]

Scale coordinates and unit-cell to get a new geometry with proper scaling

Parameters

scale (float) – the scale factor for the new geometry (lattice vectors, coordinates and the atomic radii are scaled).

set_nsc(*args, **kwargs)

Set the number of super-cells in the SuperCell object

See set_nsc for allowed parameters.

See also

SuperCell.set_nsc

the underlying called method

set_sc(sc)

Overwrites the local supercell

set_supercell(sc)

Overwrites the local supercell

sort(**kwargs)[source]

Sort atoms in a nested fashion according to various criteria

There are many ways to sort a Geometry. - by Cartesian coordinates, axis - by lattice vectors, lattice - by user defined vectors, vector - by grouping atoms, group - by a user defined function, func - by a user defined function using internal sorting algorithm, func_sort

  • a combination of the above in arbitrary order

Additionally one may sort ascending or descending.

This method allows nested sorting based on keyword arguments.

Parameters
  • atoms (list or list of list, optional) – only perform sorting algorithm for subset of atoms. This is NOT a positional dependent argument. All sorting algorithms will _only_ be performed on these atoms. If a list of indices only the given atoms will be sorted, and all other atoms will be kept in the final structure at their initial indices. If a list of list of indices, each list of indices will be sorted individually (see ret_atoms). Default, all atoms will be sorted.

  • ret_atoms (bool, optional) – return a list of list for the groups of atoms that have been sorted.

  • axis (int or tuple of int, optional) – sort coordinates according to Cartesian coordinates, if a tuple of ints is passed it will be equivalent to sort(axis0=axis[0], axis1=axis[1]). This behaves differently than numpy.lexsort!

  • lattice (int or tuple of int, optional) – sort coordinates according to lattice vectors, if a tuple of ints is passed it will be equivalent to sort(lattice0=lattice[0], lattice1=lattice[1]). Note that before sorting we multiply the fractional coordinates by the length of the lattice vector. This ensures that atol is meaningful for both axis and lattice since they will be on the same order of magnitude. This behaves differently than numpy.lexsort!

  • vector ((3, ), optional) – sort along a user defined vector, similar to lattice but with a user defined direction. Note that lattice sorting and vector sorting are only equivalent when the lattice vector is orthogonal to the other lattice vectors.

  • group ({'Z', 'symbol', 'tag', 'species'} or (str, ..), optional) – group together a set of atoms by various means. group may be one of the listed strings. For 'Z' atoms will be grouped in atomic number For 'symbol' atoms will be grouped by their atomic symbol. For 'tag' atoms will be grouped by their atomic tag. For 'species' atoms will be sorted according to their specie index. If a tuple/list is passed the first item is described. All subsequent items are a list of groups, where each group comprises elements that should be sorted on an equal footing. If one of the groups is None, that group will be replaced with all non-mentioned elements. See examples.

  • func (callable or list-like of callable, optional) – pass a sorting function which should have an interface like func(geometry, atoms, **kwargs). The first argument is the geometry to sort. The 2nd argument is a list of indices in the current group of sorted atoms. And **kwargs are any optional arguments currently collected, i.e. ascend, atol etc. The function should return either a list of atoms, or a list of list of atoms (in which case the atomic indices have been split into several groups that will be sorted individually for subsequent sorting methods). In either case the returned indices must never hold any other indices but the ones passed as atoms. If a list/tuple of functions, they will be processed in that order.

  • func_sort (callable or list-like of callable, optional) – pass a function returning a 1D array corresponding to all atoms in the geometry. The interface should simply be: func(geometry). Those values will be passed down to the internal sorting algorithm. To be compatible with atol the returned values from func_sort should be on the scale of coordinates (in Ang).

  • descend (ascend,) – control ascending or descending sorting for all subsequent sorting methods. Default ascend=True.

  • atol (float, optional) – absolute tolerance when sorting numerical arrays for subsequent sorting methods. When a selection of sorted coordinates are grouped via atol, we ensure such a group does not alter its indices. I.e. the group is always ascending indices. Note this may have unwanted side-effects if atol is very large compared to the difference between atomic coordinates. Default 1e-9.

Notes

The order of arguments is also the sorting order. sort(axis=0, lattice=0) is different from sort(lattice=0, axis=0)

All arguments may be suffixed with integers. This allows multiple keyword arguments to control sorting algorithms in different order. It also allows changing of sorting settings between different calls. Note that the integers have no relevance to the order of execution! See examples.

Returns

  • geometry (Geometry) – sorted geometry

  • index (list of list of indices) – indices that would sort the original structure to the output, only returned if ret_atoms is True

Examples

>>> geom = sisl.geom.bilayer(top_atoms=sisl.Atom[5, 7], bottom_atoms=sisl.Atom(6))
>>> geom = geom.tile(5, 0).repeat(7, 1)

Sort according to \(x\) coordinate

>>> geom.sort(axis=0)

Sort according to \(z\), then \(x\) for each group created from first sort

>>> geom.sort(axis=(2, 0))

Sort according to \(z\), then first lattice vector

>>> geom.sort(axis=2, lattice=0)

Sort according to \(z\) (ascending), then first lattice vector (descending)

>>> geom.sort(axis=2, ascend=False, lattice=0)

Sort according to \(z\) (descending), then first lattice vector (ascending) Note how integer suffixes has no importance.

>>> geom.sort(ascend1=False, axis=2, ascend0=True, lattice=0)

Sort only atoms range(1, 5) first by \(z\), then by first lattice vector

>>> geom.sort(axis=2, lattice=0, atoms=np.arange(1, 5))

Sort two groups of atoms [range(1, 5), range(5, 10)] (individually) by \(z\)

>>> geom.sort(axis=2, atoms=[np.arange(1, 5), np.arange(5, 10)])

The returned sorting indices may be used for manual sorting. Note however, that this requires one to perform a sorting for all atoms. In such a case the following sortings are equal.

>>> geom0, atoms0 = geom.sort(axis=2, lattice=0, ret_atoms=True)
>>> _, atoms1 = geom.sort(axis=2, ret_atoms=True)
>>> geom1, atoms1 = geom.sort(lattice=0, atoms=atoms1, ret_atoms=True)
>>> geom2 = geom.sub(np.concatenate(atoms0))
>>> geom3 = geom.sub(np.concatenate(atoms1))
>>> assert geom0 == geom1
>>> assert geom0 == geom2
>>> assert geom0 == geom3

Default sorting is equivalent to axis=(0, 1, 2)

>>> assert geom.sort() == geom.sort(axis=(0, 1, 2))

Sort along a user defined vector [2.2, 1., 0.]

>>> geom.sort(vector=[2.2, 1., 0.])

Integer specification has no influence on the order of operations. It is _always_ the keyword argument order that determines the operation.

>>> assert geom.sort(axis2=1, axis0=0, axis1=2) == geom.sort(axis=(1, 0, 2))

Sort by atomic numbers

>>> geom.sort(group='Z') # 5, 6, 7

One may group several elements together on an equal footing (None means all non-mentioned elements) The order of the groups are important (the first two are _not_ equal, the last three _are_ equal)

>>> geom.sort(group=('symbol', 'C', None)) # C, [B, N]
>>> geom.sort(group=('symbol', None, 'C')) # [B, N], C
>>> geom.sort(group=('symbol', ['N', 'B'], 'C')) # [B, N], C (B and N unaltered order)
>>> geom.sort(group=('symbol', ['B', 'N'], 'C')) # [B, N], C (B and N unaltered order)

A too high atol may have unexpected side-effects. This is because of the way the sorting algorithm splits the sections for nested sorting. So for coordinates with a continuous displacement the sorting may break and group a large range into 1 group. Consider the following array to be split in groups while sorting.

An example would be a linear chain with a middle point with a much shorter distance.

>>> x = np.arange(5) * 0.1
>>> x[3:] -= 0.095
y = z = np.zeros(5)
geom = si.Geometry(np.stack((x, y, z), axis=1))
>>> geom.xyz[:, 0]
[0.    0.1   0.2   0.205 0.305]

In this case a high tolerance (``atol>0.005`) would group atoms 2 and 3 together

>>> geom.sort(atol=0.01, axis=0, ret_atoms=True)[1]
[[0], [1], [2, 3], [4]]

However, a very low tolerance will not find these two as atoms close to each other.

>>> geom.sort(atol=0.001, axis=0, ret_atoms=True)[1]
[[0], [1], [2], [3], [4]]
sparserij(dtype=<class 'numpy.float64'>, na_iR=1000, method='rand')[source]

Return the sparse matrix with all distances in the matrix

The sparse matrix will only be defined for the elements which have orbitals overlapping with other atoms.

Parameters
  • dtype (numpy.dtype, numpy.float64) – the data-type of the sparse matrix

  • na_iR (int, 1000) – number of atoms within the sphere for speeding up the iter_block loop.

  • method (str, optional) – see iter_block for details

Returns

sparse matrix with all rij elements

Return type

SparseAtom

See also

iter_block

the method for looping the atoms

distance

create a list of distances

sub(atoms, cell=None)[source]

Create a new Geometry with a subset of this Geometry

Indices passed MUST be unique.

Negative indices are wrapped and thus works.

Parameters
  • atoms (int or array_like) – indices/boolean of all atoms to be removed

  • cell (array_like or SuperCell, optional) – the new associated cell of the geometry (defaults to the same cell)

See also

SuperCell.fit

update the supercell according to a reference supercell

remove

the negative of this routine, i.e. remove a subset of atoms

swap(atoms_a, atoms_b)[source]

Swap a set of atoms in the geometry and return a new one

This can be used to reorder elements of a geometry.

Parameters
  • atoms_a (array_like) – the first list of atomic coordinates

  • atoms_b (array_like) – the second list of atomic coordinates

swapaxes(axis_a, axis_b, swap='cell+xyz')[source]

Swap the axis for the atomic coordinates and the cell vectors

If swapaxes(0,1) it returns the 0 and 1 values swapped in the cell variable.

Parameters
  • axis_a (int) – axis 1, swaps with b

  • axis_b (int) – axis 2, swaps with a

  • swap ({'cell+xyz', 'cell', 'xyz'}) – decide what to swap, if ‘cell’ is in swap then the cell axis are swapped. if ‘xyz’ is in swap then the xyz (Cartesian) axis are swapped. Both may be in swap.

tile(reps, axis)[source]

Tile the geometry to create a bigger one

The atomic indices are retained for the base structure.

This method allows to utilise Bloch’s theorem when creating Hamiltonian parameter sets for TBtrans.

Tiling and repeating a geometry will result in the same geometry. The only difference between the two is the final ordering of the atoms.

Parameters
  • reps (int) – number of tiles (repetitions)

  • axis (int) – direction of tiling, 0, 1, 2 according to the cell-direction

Examples

>>> geom = Geometry([[0, 0, 0], [0.5, 0, 0]], sc=1.)
>>> g = geom.tile(2,axis=0)
>>> print(g.xyz) 
[[0.   0.   0. ]
 [0.5  0.   0. ]
 [1.   0.   0. ]
 [1.5  0.   0. ]]
>>> g = geom.tile(2,0).tile(2,axis=1)
>>> print(g.xyz) 
[[0.   0.   0. ]
 [0.5  0.   0. ]
 [1.   0.   0. ]
 [1.5  0.   0. ]
 [0.   1.   0. ]
 [0.5  1.   0. ]
 [1.   1.   0. ]
 [1.5  1.   0. ]]

See also

repeat

equivalent but different ordering of final structure

cut

opposite method of this

toASE()[source]

Returns the geometry as an ASE Atoms object

translate(v, atoms=None, cell=False)

Translates the geometry by v

One can translate a subset of the atoms by supplying atoms.

Returns a copy of the structure translated by v.

Parameters
  • v (array_like) – the vector to displace all atomic coordinates

  • atoms (int or array_like, optional) – only displace the given atomic indices, if not specified, all atoms will be displaced

  • cell (bool, optional) – If True the supercell also gets enlarged by the vector

uc2sc(atoms, unique=False)[source]

Returns atom from unit-cell indices to supercell indices, possibly removing dublicates

Parameters
  • atoms (array_like or int) – the atomic unit-cell indices to be converted to supercell indices

  • unique (bool, optional) – If True the returned indices are unique and sorted.

property volume

Returns the inherent SuperCell objects vol

within(shapes, atoms=None, atoms_xyz=None, ret_xyz=False, ret_rij=False)[source]

Indices of atoms in the entire supercell within a given shape from a given coordinate

This heavily relies on the within_sc method.

Note that if a connection is made in a neighbouring super-cell then the atomic index is shifted by the super-cell index times number of atoms. This allows one to decipher super-cell atoms from unit-cell atoms.

Parameters
  • shapes (Shape, list of Shape) –

  • atoms (array_like, optional) – List of indices for atoms that are to be considered

  • atoms_xyz (array_like, optional) – The atomic coordinates of the equivalent atoms variable (atoms must also be passed)

  • ret_xyz (bool, optional) – If true this method will return the coordinates for each of the couplings.

  • ret_rij (bool, optional) – If true this method will return the distances from the xyz_ia for each of the couplings.

Returns

  • index – indices of atoms (in supercell indices) within the shape

  • xyz – atomic coordinates of the indexed atoms (only for true ret_xyz)

  • rij – distance of the indexed atoms to the center of the shape (only for true ret_rij)

within_inf(sc, periodic=None, tol=1e-05, origo=None)[source]

Find all atoms within a provided supercell

Note this function is rather different from close and within. Specifically this routine is returning all indices for the infinite periodic system (where self.nsc > 1 or periodic is true).

Atomic coordinates lying on the boundary of the supercell will be duplicated on the neighbouring supercell images. Thus performing geom.within_inf(geom.sc) may result in more atoms than in the structure.

Notes

The name of this function may change. Currently it should only be used internally in sisl.

Parameters
  • sc (SuperCell or SuperCellChild) – the supercell in which this geometry should be expanded into.

  • periodic (list of bool) – explicitly define the periodic directions, by default the periodic directions are only where self.nsc > 1.

  • tol (float, optional) – length tolerance for the fractional coordinates to be on a duplicate site (in Ang). This allows atoms within tol of the cell boundaries to be taken as inside the cell.

  • origo ((3, ) of float) – origo that is the basis for comparison

Returns

  • numpy.ndarray – unit-cell atomic indices which are inside the sc cell

  • numpy.ndarray – atomic coordinates for the ia atoms (including supercell offsets)

  • numpy.ndarray – integer supercell offsets for ia atoms

within_sc(shapes, isc=None, atoms=None, atoms_xyz=None, ret_xyz=False, ret_rij=False)[source]

Indices of atoms in a given supercell within a given shape from a given coordinate

This returns a set of atomic indices which are within a sphere of radius R.

If R is a tuple/list/array it will return the indices: in the ranges:

>>> ( x <= R[0] , R[0] < x <= R[1], R[1] < x <= R[2] )
Parameters
  • shapes (Shape or list of Shape) –

    A list of increasing shapes that define the extend of the geometric volume that is searched. It is vital that:

    shapes[0] in shapes[1] in shapes[2] ...
    

  • isc (array_like, optional) – The super-cell which the coordinates are checked in. Defaults to [0, 0, 0]

  • atoms (array_like, optional) – List of atoms that will be considered. This can be used to only take out a certain atoms.

  • atoms_xyz (array_like, optional) – The atomic coordinates of the equivalent idx variable (idx must also be passed)

  • ret_xyz (bool, optional) – If True this method will return the coordinates for each of the couplings.

  • ret_rij (bool, optional) – If True this method will return the distance to the center of the shapes

Returns

  • index – indices of atoms (in supercell indices) within the shape

  • xyz – atomic coordinates of the indexed atoms (only for true ret_xyz)

  • rij – distance of the indexed atoms to the center of the shape (only for true ret_rij)

write(sile, *args, **kwargs)[source]

Writes geometry to the Sile using sile.write_geometry

Parameters
  • sile (Sile, str or pathlib.Path) – a Sile object which will be used to write the geometry if it is a string it will create a new sile using get_sile

  • **kwargs (*args,) –

    Any other args will be passed directly to the underlying routine

See also

read

reads a Geometry from a given Sile/file