sisl.geometry module

Geometry class to retain the atomic structure.

A Geometry contains all necessary components regarding an atomic configuration:

  1. Number of atoms
  2. Atomic coordinates (in Cartesian coordinates)
  3. Atomic species
  4. Unit cell where the atoms are contained

The class implements a wide variety of routines for manipulation of the above listed items.

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

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.

Parameters:

xyz : array_like

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

atom : array_like or Atoms

atomic species retrieved from the PeriodicTable

sc : SuperCell

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

See also

Atoms
contained atoms self.atom
Atom
contained atoms are each an object of this

Examples

An atomic lattice consisting of Hydrogen atoms. An atomic square 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
atom Atoms for the geometry (Atoms object)
maxR([all]) Maximum orbital range of the atoms
xyz (ndarray) atomic coordinates
sc (SuperCell) the supercell describing the periodicity of the geometry
no: int total number of orbitals in the geometry

Methods

ArgumentParser([parser])
a2isc(ia) Returns the super-cell index for a specific/list atom
a2o(ia[, all]) Returns an orbital index of the first orbital of said atom.
a2sc(a) Returns the super-cell offset for a specific atom
add(other) Adds atoms (as is) from the other geometry.
add_vacuum(vacuum, axis) Add vacuum along the axis lattice vector
append(other, axis) Appends structure along axis.
asc2uc(atom[, uniq]) Returns atom from super-cell indices to unit-cell indices, possibly removing dublicates
attach(s_idx, other, o_idx[, dist, axis]) Attaches another Geometry at the s_idx index with respect to o_idx using different methods.
axyz([atom, isc]) Return the atomic coordinates in the supercell of a given atom.
bond_correct(ia, atom[, method]) Corrects the bond between ia and the atom.
center([atom, which]) Returns the center of the geometry
close(xyz_ia[, R, idx, idx_xyz, ret_xyz, ...]) Returns supercell atomic indices for all atoms connecting to xyz_ia
close_all(xyz_ia[, R, idx, idx_xyz, ...]) Returns supercell atomic indices for all atoms connecting to xyz_ia
close_sc(xyz_ia[, isc, R, idx, idx_xyz, ...]) Calculates which atoms are close to some atom or point in space, only returns so relative to a super-cell.
copy() A copy of the object.
cut(seps, axis[, seg, rtol, atol]) Returns a subset of atoms from the geometry by cutting the geometry into seps parts along the direction axis.
distance([atom, R, tol, method]) Calculate the distances for all atoms in shells of radius tol within max_R
equal(other[, R]) 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, geom) 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, atom, method]) Iterator for performance critical loops
iter_block_rand([iR, R, atom]) Perform the random block-iteration by randomly selecting the next center of block
iter_block_shape([shape, iR, atom]) Perform the grid block-iteration by looping a grid
iter_orbitals([atom, local]) Returns an iterator over all atoms and their associated orbitals
iter_species([atom]) Iterator over all atoms and species as a tuple in this geometry
maxR([all]) Maximum orbital range of the atoms
mirror(plane[, atom]) Mirrors the structure around the center of the atoms
move(v[, atom, cell]) Translates the geometry by v
o2a(io) Returns an atomic index corresponding to the orbital indicies.
o2isc(io) Returns the super-cell index for a specific orbital.
o2sc(o) Returns the super-cell offset for a specific orbital.
optimize_nsc([axis, R]) Optimize the number of supercell connections based on self.maxR()
orij(io, jo) Return distance between orbital io and jo, orbitals are expected to be in super-cell indices
osc2uc(orbs[, uniq]) Returns orbitals from super-cell indices to unit-cell indices, possibly removing dublicates
prepend(other, axis) Prepends structure along axis.
read(sile, *args, **kwargs) Reads geometry from the Sile using Sile.read_geometry
remove(atom) Remove atoms from the geometry.
repeat(reps, axis) Create a repeated geometry
reverse([atom]) Returns a reversed geometry
rij(ia, ja) Distance between atom ia and ja, atoms are expected to be in super-cell indices
rotate(angle, v[, origo, atom, only, radians]) Rotate geometry around vector and return a new geometry
rotate_miller(m, v) Align Miller direction along v
rotatea(angle[, origo, atom, only, radians]) Rotate around first lattice vector
rotateb(angle[, origo, atom, only, radians]) Rotate around second lattice vector
rotatec(angle[, origo, atom, only, radians]) Rotate around third lattice vector
sc2uc(atom[, uniq]) Returns atom from super-cell 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(nsc) Set the number of super-cells in the SuperCell object
set_sc(sc) Overwrites the local supercell
set_supercell(sc) Overwrites the local supercell
sparserij([dtype, na_iR, method]) Return the sparse matrix with all distances in the matrix
sub(atom[, cell]) Create a new Geometry with a subset of this Geometry
swap(a, b) Swap a set of atoms in the geometry and return a new one
swapaxes(a, 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[, atom, cell]) Translates the geometry by v
within(shapes[, idx, idx_xyz, ret_xyz, ret_rij]) Returns supercell atomic indices for all atoms connecting to xyz_ia
within_sc(shapes[, isc, idx, idx_xyz, ...]) Calculates which atoms are close to some atom or point in space, only returns so relative to a super-cell.
write(sile, *args, **kwargs) Writes geometry to the Sile using sile.write_geometry
ArgumentParser(parser=None, *args, **kwargs)
a2isc(ia)[source]

Returns the super-cell index for a specific/list atom

Returns a vector of 3 numbers with integers.

a2o(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:

ia : 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(a)[source]

Returns the super-cell offset for a specific atom

add(other)[source]

Adds atoms (as is) from the other geometry. This will not alter the cell vectors.

Parameters:

other : Geometry

Other geometry class which is added

See also

append
appending geometries
prepend
prending geometries
attach
attach a geometry
insert
insert a geometry
append(other, axis)[source]

Appends structure 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.

See also

add
add geometries
prepend
prending geometries
attach
attach a geometry
insert
insert a geometry
asc2uc(atom, uniq=False)

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

atom

Atoms for the geometry (Atoms object)

atoms

Atoms for the geometry (Atoms object)

attach(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.

Parameters:

dist : array_like, float, str ('calc')

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.

axis : int

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

axyz(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:

atom : 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(cell=1., xyz=[[0,0,0],[0.5,0,0]])
>>> print(geom.axyz(isc=[1,0,0])
[[ 1.   0.   0. ]
 [ 1.5  0.   0. ]]
>>> geom = Geometry(cell=1., xyz=[[0,0,0],[0.5,0,0]])
>>> print(geom.axyz(0))
[ 1.   0.   0. ]
bond_correct(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:

ia : int

The atom to be displaced according to the atomic radius

atom : array_like or int

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

method : str, float

If str will use that as lookup in Atom.radius. Else it will be the new bond-length.

center(atom=None, which='xyz')[source]

Returns the center of the geometry

By specifying which one can control whether it should be:

  • xyz|position: Center of coordinates (default)
  • mass: Center of mass
  • cell: Center of cell
Parameters:

atom : array_like

list of atomic indices to find center of

which : {‘xyz’, ‘mass’, ‘cell’}

determine whether center should be of ‘cell’, mass-centered (‘mass’), or absolute center of the positions.

close(xyz_ia, R=None, idx=None, idx_xyz=None, ret_xyz=False, ret_rij=False)[source]

Returns supercell atomic indices for all atoms connecting to xyz_ia

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``

idx : array_like, optional

List of indices for atoms that are to be considered

idx_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 distances from the xyz_ia for each of the couplings.

close_all(xyz_ia, R=None, idx=None, idx_xyz=None, ret_xyz=False, ret_rij=False)

Returns supercell atomic indices for all atoms connecting to xyz_ia

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``

idx : array_like, optional

List of indices for atoms that are to be considered

idx_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 distances from the xyz_ia for each of the couplings.

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

Calculates which atoms are close to some atom or point in space, only returns so relative to a super-cell.

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

idx : array_like of int, optional

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

idx_xyz : array_like of float, 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 for each of the couplings.

copy()[source]

A copy of the object.

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

Returns a subset of atoms from the geometry by cutting the geometry into seps parts along the direction axis. It will then _only_ return the first cut.

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.

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)

distance(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:

atom : 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.

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()
[ 1.]
>>> geom.distance(tol=[0.5, 0.4, 0.3, 0.2])
[ 1.          1.41421356]
>>> geom.distance(R=2, tol=[0.5, 0.4, 0.3, 0.2])
[ 1.          1.41421356  2.        ]
>>> geom.distance(R=2, tol=[0.5, 0.7]) # the R = 1 and R = 2 ** .5 gets averaged
[ 1.20710678  2.        ]
equal(other, R=True)[source]

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

Parameters:

other : Geometry

the other Geometry to check against

maxR : bool, optional

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

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 sisl package object will be created.

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

insert(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:

atom : int

the index at which atom the other geometry is inserted

geom : Geometry

the other geometry to be inserted

See also

add
add geometries
append
appending geometries
prepend
prending geometries
attach
attach a geometry
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=10, 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:

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

atom : 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:

Two lists with [0] being a list of atoms to be looped and [1] being the atoms that

need searched.

iter_block_rand(iR=10, R=None, atom=None)[source]

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

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

Perform the grid block-iteration by looping a grid

iter_orbitals(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:

atom : 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.

See also

iter
iterate over atomic indices
iter_species
iterate across indices and atomic species
iter_species(atom=None)[source]

Iterator over all atoms 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:

atom : 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
lasto

The last orbital on the corresponding atom

mass

Returns the mass of all atoms as an array

maxR(all=False)[source]

Maximum orbital range of the atoms

mirror(plane, atom=None)[source]

Mirrors the structure around the center of the atoms

move(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:

v : array_like

the vector to displace all atomic coordinates

atom : 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

na

Number of atoms in geometry

na_s

Number of supercell atoms

no

Number of orbitals

no_s

Number of supercell orbitals

o2a(io)[source]

Returns an 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

o2isc(io)[source]

Returns the super-cell index for a specific orbital.

Returns a vector of 3 numbers with integers.

o2sc(o)[source]

Returns the super-cell offset for a specific orbital.

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

orbitals

List of orbitals per atom

orij(io, jo)[source]

Return distance between orbital io and jo, orbitals are expected to be in super-cell indices

Returns the distance between two orbitals:

\[\begin{split}r\\_{ij} = |r\\_j - r\\_i|\end{split}\]
Parameters:

io : int or array_like

orbital index of first orbital

jo : int or array_like

orbital indices

osc2uc(orbs, uniq=False)[source]

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

prepend(other, axis)[source]

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

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

See also

add
add geometries
append
appending geometries
attach
attach a geometry
insert
insert a geometry
static read(sile, *args, **kwargs)[source]

Reads geometry from the Sile using Sile.read_geometry

Parameters:

sile : Sile, str

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.

remove(atom)[source]

Remove atoms from the geometry.

Indices passed MUST be unique.

Negative indices are wrapped and thus works.

Parameters:

atom : array_like

indices of all atoms to be removed.

repeat(reps, axis)[source]

Create a repeated geometry

The atomic indices are NOT retained for 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.

It is adviced to only use this for electrode Bloch’s theorem purposes as tile is faster.

Parameters:

reps : int

number of repetitions

axis : int

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

See also

tile
equivalent but different ordering of final structure

Examples

>>> geom = Geometry(cell=[[1.,0,0],[0,1.,0.],[0,0,1.]],xyz=[[0,0,0],[0.5,0,0]])
>>> 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. ]
 [ 1.   0.   0. ]
 [ 0.   1.   0. ]
 [ 1.   1.   0. ]
 [ 0.5  0.   0. ]
 [ 1.5  0.   0. ]
 [ 0.5  1.   0. ]
 [ 1.5  1.   0. ]]
reverse(atom=None)[source]

Returns a reversed geometry

Also enables reversing a subset of the atoms.

Parameters:

atom : 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 are expected to be in super-cell indices

Returns the distance between two atoms:

\[\begin{split}r\\_{ij} = |r\\_j - r\\_i|\end{split}\]
Parameters:

ia : int or array_like

atomic index of first atom

ja : int or array_like

atomic indices

rotate(angle, v, origo=None, atom=None, only='abc+xyz', radians=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 radians of which the geometry should be rotated

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.

atom : 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

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.

rotatea(angle, origo=None, atom=None, only='abc+xyz', radians=False)[source]

Rotate around first lattice vector

See also

rotate
called routine with v = self.cell[0, :]
rotateb(angle, origo=None, atom=None, only='abc+xyz', radians=False)[source]

Rotate around second lattice vector

See also

rotate
called routine with v = self.cell[1, :]
rotatec(angle, origo=None, atom=None, only='abc+xyz', radians=False)[source]

Rotate around third lattice vector

See also

rotate
called routine with v = self.cell[2, :]
sc2uc(atom, uniq=False)[source]

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

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

sparserij(dtype=<type '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, ‘rand’

see iter_block for details

Returns:

SparseCSR

sparse matrix with all rij elements

See also

iter_block
the method for looping the atoms
distance
create a list of distances
sub(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:

atom : array_like

indices 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
swap(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:

a : array_like

the first list of atomic coordinates

b : array_like

the second list of atomic coordinates

swapaxes(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:

a : int

axes 1, swaps with b

b : int

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(reps, axis)[source]

Tile the geometry to create a bigger one

The atomic indices are retained for the base structure.

Parameters:

reps : int

number of tiles (repetitions)

axis : int

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

See also

repeat
equivalent but different ordering of final structure

Examples

>>> geom = Geometry(cell=[[1.,0,0],[0,1.,0.],[0,0,1.]],xyz=[[0,0,0],[0.5,0,0]])
>>> 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. ]]
toASE()[source]

Returns the geometry as an ASE Atoms object

translate(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:

v : array_like

the vector to displace all atomic coordinates

atom : 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

within(shapes, idx=None, idx_xyz=None, ret_xyz=False, ret_rij=False)[source]

Returns supercell atomic indices for all atoms connecting to xyz_ia

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:

shapes : Shape, list of Shape

idx : array_like, optional

List of indices for atoms that are to be considered

idx_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 distances from the xyz_ia for each of the couplings.

within_sc(shapes, isc=None, idx=None, idx_xyz=None, ret_xyz=False, ret_rij=False)[source]

Calculates which atoms are close to some atom or point in space, only returns so relative to a super-cell.

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]

idx : array_like, optional

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

idx_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 for each of the couplings.

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

Writes geometry to the Sile using sile.write_geometry

Parameters:

sile : Sile, str

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

sisl.geometry.sgeom(geom=None, argv=None, ret_geometry=False)[source]

Main script for sgeom script.

This routine may be called with argv and/or a Sile which is the geometry at hand.

Parameters:

geom : Geometry, BaseSile

this may either be the geometry, as-is, or a Sile which contains the geometry.

argv : list of str

the arguments passed to sgeom

ret_geometry : bool (False)

whether the function should return the geometry