Geometry

class sisl.Geometry(xyz, atom=None, sc=None, names=None)[source]

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]}
}
Parameters
xyzarray_like

atomic coordinates xyz[i, :] is the atomic coordinate of the i’th atom.

atomarray_like or Atoms

atomic species retrieved from the PeriodicTable

scSuperCell

the unit-cell describing the atoms in a periodic super-cell

See also

Atoms

contained atoms self.atoms

Atom

contained atoms are each an object of this

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'))
Attributes
na

Number of atoms in geometry

xyzndarray

atomic coordinates

atoms

Atoms for the geometry (Atoms object)

orbitals

List of orbitals per atom

scSuperCell

the supercell describing the periodicity of the geometry

no

Number of orbitals

n_sint

Returns the inherent SuperCell objects n_s

no_sint

Number of supercell orbitals

Attributes

atom

Atoms for the geometry (Atoms object)

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(self, ia, ja)

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

__init__(self, xyz[, atom, sc, names])

Initialize self.

a2isc(self, ia)

Returns super-cell index for a specific/list atom

a2o(self, ia[, all])

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

a2sc(self, a)

Returns the super-cell offset for a specific atom

a2transpose(self, atom1[, atom2])

Transposes connections from atom1 to atom2 such that supercell connections are transposed

add(self, other)

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

add_vacuum(self, vacuum, axis)

Add vacuum along the axis lattice vector

angle(self, atom[, dir, ref, rad])

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

append(self, other, axis[, align])

Appends two structures along axis

area(self, ax0, ax1)

Calculate the area spanned by the two axis ax0 and ax1

as_primary(self, na_primary[, ret_super])

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

asc2uc(self, atom[, unique])

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

attach(self, s_idx, other, o_idx[, dist, axis])

Attaches another Geometry at the s_idx index with respect to o_idx using different methods.

auc2sc(self, atom[, unique])

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

axyz(self[, atom, isc])

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

bond_correct(self, ia, atom[, method])

Corrects the bond between ia and the atom.

center(self[, atom, what])

Returns the center of the geometry

close(self, xyz_ia[, R, idx, idx_xyz, …])

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

close_sc(self, xyz_ia[, isc, R, idx, …])

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

copy(self)

A copy of the object.

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

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

distance(self[, atom, R, tol, method])

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

equal(self, 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(self[, na, iR, R])

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

insert(self, atom, geom)

Inserts other atoms right before index

is_orthogonal(self)

Return true if all cell vectors are linearly independent

iter(self)

An iterator over all atomic indices

iter_block(self[, iR, R, atom, method])

Iterator for performance critical loops

iter_block_rand(self[, iR, R, atom])

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

iter_block_shape(self[, shape, iR, atom])

Perform the grid block-iteration by looping a grid

iter_orbitals(self[, atom, local])

Returns an iterator over all atoms and their associated orbitals

iter_species(self[, atom])

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

maxR(self[, all])

Maximum orbital range of the atoms

mirror(self, plane[, atom])

Mirrors the atomic coordinates by multiplying by -1

move(self, v[, atom, cell])

Translates the geometry by v

o2a(self, io[, unique])

Atomic index corresponding to the orbital indicies.

o2isc(self, io)

Returns the super-cell index for a specific orbital.

o2sc(self, o)

Returns the super-cell offset for a specific orbital.

o2transpose(self, orb1[, orb2])

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

oRij(self, io, jo)

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

optimize_nsc(self[, axis, R])

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

orij(self, io, jo)

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

osc2uc(self, orb[, unique])

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

ouc2sc(self, orb[, unique])

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

prepend(self, other, axis[, align])

Prepend two structures along axis

read(sile, \*args, \*\*kwargs)

Reads geometry from the Sile using Sile.read_geometry

reduce(self)

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

remove(self, atom)

Remove atoms from the geometry.

reorder(self)

Reorders atoms according to first occurence in the geometry

repeat(self, reps, axis)

Create a repeated geometry

reverse(self[, atom])

Returns a reversed geometry

rij(self, ia, ja)

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

rotate(self, angle, v[, origo, atom, only, rad])

Rotate geometry around vector and return a new geometry

rotate_miller(self, m, v)

Align Miller direction along v

sc2uc(self, atom[, unique])

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

sc_index(self, \*args, \*\*kwargs)

Call local SuperCell object sc_index function

scale(self, scale)

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

set_nsc(self, \*args, \*\*kwargs)

Set the number of super-cells in the SuperCell object

set_sc(self, sc)

Overwrites the local supercell

set_supercell(self, sc)

Overwrites the local supercell

sort(self[, axes])

Return an equivalent geometry by sorting the coordinates according to the axis orders

sparserij(self[, dtype, na_iR, method])

Return the sparse matrix with all distances in the matrix

sub(self, atom[, cell])

Create a new Geometry with a subset of this Geometry

swap(self, a, b)

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

swapaxes(self, a, b[, swap])

Swap the axis for the atomic coordinates and the cell vectors

tile(self, reps, axis)

Tile the geometry to create a bigger one

toASE(self)

Returns the geometry as an ASE Atoms object

translate(self, v[, atom, cell])

Translates the geometry by v

uc2sc(self, atom[, unique])

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

within(self, shapes[, idx, idx_xyz, …])

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

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

Find all atoms within a provided supercell

within_sc(self, shapes[, isc, idx, idx_xyz, …])

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

write(self, sile, \*args, \*\*kwargs)

Writes geometry to the Sile using sile.write_geometry

Rij(self, 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
iaint or array_like

atomic index of first atom

jaint or array_like

atomic indices

a2isc(self, ia)[source]

Returns super-cell index for a specific/list atom

Returns a vector of 3 numbers with integers.

a2o(self, ia, 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
iaarray_like

Atomic indices

allbool, optional

False, return only the first orbital corresponding to the atom, True, returns list of the full atom

a2sc(self, a)[source]

Returns the super-cell offset for a specific atom

a2transpose(self, atom1, atom2=None)[source]

Transposes connections from atom1 to atom2 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.

Parameters
atom1array_like

atomic indices must have same length as atom2 or length 1

atom2array_like, optional

atomic indices must have same length as atom1 or length 1. If not present then only atom1 will be returned in transposed indices.

Returns
atom2array_like

transposed indices for atom2 (only returned if atom2 is not None)

atom1array_like

transposed indices for atom1

Examples

>>> gr = geom.graphene()
>>> idx = gr.close(0, 1.5)
>>> idx
array([0, 1, 5, 9], dtype=int32)
>>> gr.a2transpose(0, idx)
(array([0, 1, 1, 1], dtype=int32), array([ 0,  0, 14, 10], dtype=int32))
add(self, other)[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
otherGeometry or SuperCell

Other geometry class which is added

See also

append

appending geometries

prepend

prending geometries

attach

attach a geometry

insert

insert a geometry

add_vacuum(self, vacuum, axis)

Add vacuum along the axis lattice vector

Parameters
vacuumfloat

amount of vacuum added, in Ang

axisint

the lattice vector to add vacuum along

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

The angle between atom atom 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
atomint or array_like

indices/boolean of all atoms to be removed

dirstr, int or vector

the direction from which the angle is calculated from, default to x

refint or coordinate, optional

the reference point from which the vectors are drawn, default to origo

radbool, optional

whether the returned value is in radians

append(self, other, axis, align='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
otherGeometry or SuperCell

Other geometry class which needs to be appended If a SuperCell only the super cell will be extended

axisint

Cell direction to which the other geometry should be appended.

align{‘none’, ‘min’}

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.

See also

add

add geometries

prepend

prending geometries

attach

attach a geometry

insert

insert a geometry

area(self, ax0, ax1)

Calculate the area spanned by the two axis ax0 and ax1

as_primary(self, 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_primaryint

number of atoms in the primary unit cell

ret_superbool, 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
SislErrorin case the algorithm fails.
asc2uc(self, atom, unique=False)

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

Parameters
atomarray_like or int

the atomic supercell indices to be converted to unit-cell indices

uniquebool, optional

If True the returned indices are unique and sorted.

property atom

Atoms for the geometry (Atoms object)

property atoms

Atoms for the geometry (Atoms object)

attach(self, s_idx, other, o_idx, dist='calc', axis=None)[source]

Attaches another Geometry at the s_idx index with respect to o_idx using different methods.

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

Parameters
s_idxint

atomic index which is the base position of the attachment. The distance between s_idx and o_idx is dist.

otherGeometry

the other Geometry to attach at the given point. In this case dist from s_idx.

o_idxint

the index of the atom in other that is inserted at s_idx.

distarray_like or float or str, optional

the distance (in Ang) between the attached coordinates. If dist is arraylike 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.

axisint

specify the direction of the lattice vectors used. Not used if dist is an array-like argument.

auc2sc(self, atom, unique=False)

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

Parameters
atomarray_like or int

the atomic unit-cell indices to be converted to supercell indices

uniquebool, optional

If True the returned indices are unique and sorted.

axyz(self, atom=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
atomint or array_like

atom(s) from which we should return the coordinates, the atomic indices may be in supercell format.

iscarray_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])) # doctest: +NORMALIZE_WHITESPACE
[[1.   0.   0. ]
 [1.5  0.   0. ]]
>>> geom = Geometry([[0, 0, 0], [0.5, 0, 0]], sc=1.)
>>> print(geom.axyz(0)) # doctest: +NORMALIZE_WHITESPACE
[0.  0.  0.]
bond_correct(self, ia, atom, method='calc')[source]

Corrects the bond between ia and the atom.

Corrects the bond-length between atom ia and atom 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
iaint

The atom to be displaced according to the atomic radius

atomarray_like or int

The atom(s) from which the radius should be reduced.

methodstr, 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(self, atom=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
atomarray_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(self, xyz_ia, R=None, idx=None, idx_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_iacoordinate/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
idxarray_like, optional

List of indices for atoms that are to be considered

idx_xyzarray_like, optional

The atomic coordinates of the equivalent idx variable (idx must also be passed)

ret_xyzbool, optional

If true this method will return the coordinates for each of the couplings.

ret_rijbool, 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(self, xyz_ia, isc=(0, 0, 0), R=None, idx=None, idx_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_iaarray_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,:]).

iscarray_like, optional

The super-cell which the coordinates are checked in.

Rfloat 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.

idxarray_like of int, optional

List of atoms that will be considered. This can be used to only take out a certain atoms.

idx_xyzarray_like of float, optional

The atomic coordinates of the equivalent idx variable (idx must also be passed)

ret_xyzbool, optional

If True this method will return the coordinates for each of the couplings.

ret_rijbool, 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(self)[source]

A copy of the object.

cut(self, 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
sepsint

number of times the structure will be cut.

axisint

the axis that will be cut

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

See also

tile

opposite method of this

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
distance(self, atom=None, R=None, tol=0.1, method='average')[source]

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

Parameters
atomint or array_like, optional

only create list of distances from the given atoms, default to all atoms

Rfloat, 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.

tolfloat 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.

Returns
numpy.ndarray

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

See also

sparserij

return a sparse matrix will all distances between atoms

Examples

>>> geom = Geometry([0]*3, Atom(1, R=1.), sc=SuperCell(1., nsc=[5, 5, 1]))
>>> geom.distance() # use geom.maxR() # doctest: +NORMALIZE_WHITESPACE
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]) # doctest: +NORMALIZE_WHITESPACE
array([1.        ,  1.41421356,  2.        ])
>>> geom.distance(R=2, tol=[0.5, 0.7]) # the R = 1 and R = 2 ** .5 gets averaged # doctest: +NORMALIZE_WHITESPACE
array([1.20710678,  2.        ])
equal(self, other, R=True, tol=0.0001)[source]

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

Parameters
otherGeometry

the other Geometry to check against

Rbool, optional

if True also check if the orbital radii are the same (see Atom.equal)

tolfloat, 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
asegASE 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(self, na=1000, iR=20, R=None)[source]

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

Parameters
naint, optional

number of atoms within the radius

iRint, optional

initial iR value, which the sphere is estitametd from

Rfloat, optional

the value used for atomic range (defaults to self.maxR())

Returns
int

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

property icell

Returns the inherent SuperCell objects icell

insert(self, atom, geom)[source]

Inserts other atoms right before index

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

Parameters
atomint

the index at which atom the other geometry is inserted

geomGeometry

the other geometry to be inserted

See also

add

add geometries

append

appending geometries

prepend

prending geometries

attach

attach a geometry

is_orthogonal(self)

Return true if all cell vectors are linearly independent

property isc_off

Returns the inherent SuperCell objects isc_off

iter(self)[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(self, iR=20, R=None, atom=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
iRint, optional

the number of R ranges taken into account when doing the iterator

Rfloat, optional

enables overwriting the local R quantity. Defaults to self.maxR()

atomarray_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(self, iR=20, R=None, atom=None)[source]

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

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

Perform the grid block-iteration by looping a grid

iter_orbitals(self, atom=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
atomint or array_like, optional

only loop on the given atoms, default to all atoms

localbool, optional

whether the orbital index is the global index, or the local index relative to the atom it resides on.

See also

iter

iterate over atomic indices

iter_species

iterate across indices and atomic species

iter_species(self, atom=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
atomint 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(self, all=False)[source]

Maximum orbital range of the atoms

mirror(self, plane, atom=None)[source]

Mirrors the atomic coordinates by multiplying by -1

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

Parameters
plane{‘xy’/’ab’, ‘yz’/’bc’, ‘xz’/’ac’}

mirror the structure across the lattice vector plane

atomarray_like, optional

only mirror a subset of atoms

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

Translates the geometry by v

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

Returns a copy of the structure translated by v.

Parameters
varray_like

the vector to displace all atomic coordinates

atomint or array_like, optional

only displace the given atomic indices, if not specified, all atoms will be displaced

cellbool, 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(self, 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
ioarray_like

List of indices to return the atoms for

uniquebool, optional

If True only return the unique atoms.

o2isc(self, io)[source]

Returns the super-cell index for a specific orbital.

Returns a vector of 3 numbers with integers.

o2sc(self, o)[source]

Returns the super-cell offset for a specific orbital.

o2transpose(self, 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.

Parameters
orb1array_like

orbital indices must have same length as orb2 or length 1

orb2array_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
orb2array_like

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

orb1array_like

transposed indices for orb1

Examples

>>> gr = geom.graphene() # one orbital per site
>>> idx = gr.close(0, 1.5)
>>> idx
array([0, 1, 5, 9], dtype=int32)
>>> gr.o2transpose(0, idx)
(array([0, 1, 1, 1], dtype=int32), array([ 0,  0, 14, 10], dtype=int32))
oRij(self, 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
ioint or array_like

orbital index of first orbital

joint or array_like

orbital indices

optimize_nsc(self, 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
axisint or array_like, optional

only optimize the specified axis (default to all)

Rfloat, optional

the maximum connection radius for each atom

property orbitals

List of orbitals per atom

property origo

Returns the inherent SuperCell objects origo

orij(self, 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
ioint or array_like

orbital index of first orbital

joint or array_like

orbital indices

osc2uc(self, orb, unique=False)[source]

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

Parameters
orbarray_like or int

the orbital supercell indices to be converted to unit-cell indices

uniquebool, optional

If True the returned indices are unique and sorted.

ouc2sc(self, orb, unique=False)[source]

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

Parameters
orbarray_like or int

the orbital unit-cell indices to be converted to supercell indices

uniquebool, optional

If True the returned indices are unique and sorted.

prepend(self, other, axis, align='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
otherGeometry or SuperCell

Other geometry class which needs to be prepended If a SuperCell only the super cell will be extended

axisint

Cell direction to which the other geometry should be prepended

align{‘none’, ‘min’}

By default prepending 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.

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
sileSile, 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(self)[source]

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

Notes

This is an in-place operation.

remove(self, atom)[source]

Remove atoms from the geometry.

Indices passed MUST be unique.

Negative indices are wrapped and thus works.

Parameters
atomint 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(self)[source]

Reorders atoms according to first occurence in the geometry

Notes

This is an in-place operation.

repeat(self, 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
repsint

number of repetitions

axisint

direction of repetition, 0, 1, 2 according to the cell-direction

See also

tile

equivalent but different ordering of final structure

Examples

>>> geom = Geometry([[0, 0, 0], [0.5, 0, 0]], sc=1)
>>> g = geom.repeat(2,axis=0)
>>> print(g.xyz) # doctest: +NORMALIZE_WHITESPACE
[[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) # doctest: +NORMALIZE_WHITESPACE
[[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. ]]
reverse(self, atom=None)[source]

Returns a reversed geometry

Also enables reversing a subset of the atoms.

Parameters
atomint or array_like, optional

only reverse the given atomic indices, if not specified, all atoms will be reversed

rij(self, 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
iaint or array_like

atomic index of first atom

jaint or array_like

atomic indices

rotate(self, angle, v, origo=None, atom=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
anglefloat

the angle in degrees to rotate the geometry. Set the rad argument to use radians.

varray_like

the normal vector to the rotated plane, i.e. v = [1,0,0] will rotate the yz plane

origoint 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.

atomint 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

radbool, optional

if True the angle is provided in radians (rather than degrees)

See also

Quaternion

class to rotate

rotate_miller(self, m, v)[source]

Align Miller direction along v

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

sc2uc(self, atom, unique=False)[source]

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

Parameters
atomarray_like or int

the atomic supercell indices to be converted to unit-cell indices

uniquebool, optional

If True the returned indices are unique and sorted.

sc_index(self, *args, **kwargs)

Call local SuperCell object sc_index function

property sc_off

Returns the inherent SuperCell objects sc_off

scale(self, scale)[source]

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

Parameters
scalefloat

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

set_nsc(self, *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(self, sc)

Overwrites the local supercell

set_supercell(self, sc)

Overwrites the local supercell

sort(self, axes=(2, 1, 0))[source]

Return an equivalent geometry by sorting the coordinates according to the axis orders

Parameters
axestuple, optional

sorting axes (note the last element has highest precedence)

Returns
Geometry

sorted geometry

Examples

>>> idx = np.lexsort((self.xyz[:, i] for i in axis))
>>> new = self.sub(idx)
sparserij(self, 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
dtypenumpy.dtype, numpy.float64

the data-type of the sparse matrix

na_iRint, 1000

number of atoms within the sphere for speeding up the iter_block loop.

methodstr, optional

see iter_block for details

Returns
SparseAtom

sparse matrix with all rij elements

See also

iter_block

the method for looping the atoms

distance

create a list of distances

sub(self, atom, 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
atomint or array_like

indices/boolean of all atoms to be removed

cellarray_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(self, a, 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
aarray_like

the first list of atomic coordinates

barray_like

the second list of atomic coordinates

swapaxes(self, a, 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
aint

axes 1, swaps with b

bint

axes 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(self, 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
repsint

number of tiles (repetitions)

axisint

direction of tiling, 0, 1, 2 according to the cell-direction

See also

repeat

equivalent but different ordering of final structure

cut

opposite method of this

Examples

>>> geom = Geometry([[0, 0, 0], [0.5, 0, 0]], sc=1.)
>>> g = geom.tile(2,axis=0)
>>> print(g.xyz) # doctest: +NORMALIZE_WHITESPACE
[[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) # doctest: +NORMALIZE_WHITESPACE
[[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. ]]
toASE(self)[source]

Returns the geometry as an ASE Atoms object

translate(self, v, atom=None, cell=False)

Translates the geometry by v

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

Returns a copy of the structure translated by v.

Parameters
varray_like

the vector to displace all atomic coordinates

atomint or array_like, optional

only displace the given atomic indices, if not specified, all atoms will be displaced

cellbool, optional

If True the supercell also gets enlarged by the vector

uc2sc(self, atom, unique=False)[source]

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

Parameters
atomarray_like or int

the atomic unit-cell indices to be converted to supercell indices

uniquebool, optional

If True the returned indices are unique and sorted.

property volume

Returns the inherent SuperCell objects vol

within(self, shapes, idx=None, idx_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
shapesShape, list of Shape
idxarray_like, optional

List of indices for atoms that are to be considered

idx_xyzarray_like, optional

The atomic coordinates of the equivalent idx variable (idx must also be passed)

ret_xyzbool, optional

If true this method will return the coordinates for each of the couplings.

ret_rijbool, 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(self, 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.

Parameters
scSuperCell or SuperCellChild

the supercell in which this geometry should be expanded into.

periodiclist of bool

explicitly define the periodic directions, by default the periodic directions are only where self.nsc > 1.

tolfloat, 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

Notes

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

within_sc(self, shapes, isc=None, idx=None, idx_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
shapesShape 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] ...
iscarray_like, optional

The super-cell which the coordinates are checked in. Defaults to [0, 0, 0]

idxarray_like, optional

List of atoms that will be considered. This can be used to only take out a certain atoms.

idx_xyzarray_like, optional

The atomic coordinates of the equivalent idx variable (idx must also be passed)

ret_xyzbool, optional

If True this method will return the coordinates for each of the couplings.

ret_rijbool, 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(self, sile, *args, **kwargs)[source]

Writes geometry to the Sile using sile.write_geometry

Parameters
sileSile, 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

*args, **kwargs:

Any other args will be passed directly to the underlying routine

See also

read

reads a Geometry from a given Sile/file