lattpy.structure

Lattice structure object for defining the atom basis and neighbor connections.

class lattpy.structure.LatticeStructure(basis, **kwargs)[source]

Bases: LatticeBasis

Structure object representing a infinite Bravais lattice.

Combines the LatticeBasis with a set of atoms and connections between them to define a general lattice structure.

Inheritance

Inheritance diagram of LatticeStructure

Parameters
basisarray_like or float or LatticeBasis

The primitive basis vectors that define the unit cell of the lattice. If a LatticeBasis instance is passed it is copied and used as the new basis of the lattice.

**kwargs

Key-word arguments. Used for quickly configuring a LatticeStructure instance. Allowed keywords are:

Properties: atoms: Dictionary containing the atoms to add to the lattice. cons: Dictionary containing the connections to add to the lattice.

Examples

Two dimensional lattice with one atom in the unit cell and nearest neighbors

>>> import lattpy as lp
>>> latt = lp.LatticeStructure(np.eye(2))
>>> latt.add_atom()
>>> latt.add_connections(1)
>>> latt
LatticeStructure(dim: 2, num_base: 1, num_neighbors: [4])

Quick-setup of the same lattice:

>>> import lattpy as lp
>>> import matplotlib.pyplot as plt
>>> latt = lp.LatticeStructure.square(atoms={(0.0, 0.0): "A"}, cons={("A", "A"): 1})
>>> _ = latt.plot_cell()
>>> plt.show()

(png, pdf)

../_images/lattpy-structure-1.png
DIST_DECIMALS = 6
property num_base

int: The number of atoms in the unit cell.

property atoms

list of Atom: List of the atoms in the unit cell.

property atom_positions

list of np.ndarray: List of positions of the atoms in the unit cell.

property num_distances

int: The maximal number of distances between the lattice sites.

property num_neighbors

np.ndarray: The number of neighbors of each atom in the unitcell.

property base_neighbors

np.ndarray: The neighbors of the unitcell at the origin.

property distances

np.ndarray: List of distances between the lattice sites.

add_atom(pos=None, atom=None, primitive=False, neighbors=0, **kwargs)[source]

Adds a site to the basis of the lattice unit cell.

Parameters
pos(N) array_like or float, optional

Position of site in the unit-cell. The default is the origin of the cell. The size of the array has to match the dimension of the lattice.

atomstr or dict or Atom, optional

Identifier of the site. If a string is passed, a new Atom instance is created.

primitivebool, optional

Flag if the specified position is in cartesian or lattice coordinates. If True the passed position will be multiplied with the lattice vectors. The default is False (cartesian coordinates).

neighborsint, optional

The number of neighbor distance to calculate. If the number is 0 the distances have to be calculated manually after configuring the lattice basis.

**kwargs

Keyword arguments for ´Atom´ constructor. Only used if a new Atom instance is created.

Returns
atomAtom

The Atom-instance of the newly added site.

Raises
ValueError

Raised if the dimension of the position does not match the dimension of the lattice.

ConfigurationError

Raised if the position of the new atom is already occupied.

Examples

Construct a square lattice

>>> latt = LatticeStructure(np.eye(2))

Create an atom and add it to the origin of the unit cell of the lattice

>>> atom1 = Atom(name="A")
>>> latt.add_atom([0.0, 0.0], atom=atom1)
>>> latt.get_atom(0)
Atom(A, size=10, 0)

An Atom instance can also be created by passing the name of the (new) atom and optional keyword arguments for the constructor:

>>> latt.add_atom([0.5, 0.5], atom="B", size=15)
>>> latt.get_atom(1)
Atom(B, size=15, 1)
get_alpha(atom)[source]

Returns the index of the atom in the unit-cell.

Parameters
atomint or str or Atom

The argument for getting the atom. If a int is passed it is interpreted as the index, if a str is passed as the name of an atom.

Returns
alphaint or list of int

The atom indices. If a string was passed multiple atoms with the same name can be returned as list.

Examples

Construct a lattice with two identical atoms and a third atom in the unit cell:

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom([0, 0], atom="A")
>>> latt.add_atom([0.5, 0], atom="B")
>>> latt.add_atom([0.5, 0.5], atom="B")

Get the atom index of atom A:

>>> latt.get_alpha("A")
[0]

Get the indices of the atoms B:

>>> latt.get_alpha("B")
[1, 2]

Since there are two atoms B in the unit cell of the lattice both indices are returned.

get_atom(atom)[source]

Returns the Atom instance of the given atom in the unit cell.

Parameters
atomint or str or Atom

The argument for getting the atom. If a int is passed it is interpreted as the index, if a str is passed as the name of an atom.

Returns
atomAtom

The Atom instance of the given site.

See also

get_alpha

Get the index of the given atom.

Examples

Construct a lattice with one atom in the unit cell

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom([0, 0], atom="A")

Get the atom instance by the name

>>> latt.get_atom("A")
Atom(A, size=10, 0)

or by the index:

>>> latt.get_atom(0)
Atom(A, size=10, 0)
add_connection(atom1, atom2, num_distances=1, analyze=False)[source]

Sets the number of distances for a specific connection between two atoms.

Parameters
atom1int or str or Atom

The first atom of the connected pair.

atom2int or str or Atom

The second atom of the connected pair.

num_distancesint, optional

The number of neighbor-distance levels, e.g. setting to 1 means only nearest neighbors. The default are nearest neighbor connections.

analyzebool, optional

If True the lattice basis is analyzed after adding connections. If False the analyze-method needs to be called manually. The default is False.

See also

add_connections

Set up connections for all atoms in the unit cell in one call.

analyze

Called after setting up all the lattice connections.

Examples

Construct a square lattice with two atoms, A and B, in the unit cell:

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom([0.0, 0.0], atom="A")
>>> latt.add_atom([0.5, 0.5], atom="B")

Set next nearest and nearest neighbors between the A atoms:

>>> latt.add_connection("A", "A", num_distances=2)

Set nearest neighbors between A and B:

>>> latt.add_connection("A", "B", num_distances=1)

Set nearest neighbors between the B atoms:

>>> latt.add_connection("B", "B", num_distances=1)
add_connections(num_distances=1, analyze=True)[source]

Sets the number of distances for all possible atom-pairs of the unitcell.

Parameters
num_distancesint, optional

The number of neighbor-distance levels, e.g. setting to 1 means only nearest neighbors. The default are nearest neighbor connections.

analyzebool, optional

If True the lattice basis is analyzed after adding connections. If False the analyze-method needs to be called manually. The default is True.

See also

add_connection

Set up connection between two specific atoms in the unit cell.

analyze

Called after setting up all the lattice connections.

Examples

Construct a square lattice with one atom in the unit cell:

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom()

Set nearest neighbor hopping:

>>> latt.add_connections(num_distances=1)
analyze()[source]

Analyzes the structure of the lattice and stores neighbor data.

Check’s distances between all sites of the bravais lattice and saves the n lowest values. The neighbor lattice-indices of the unit-cell are also stored for later use. This speeds up many calculations like finding nearest neighbors.

Raises
NoAtomsError

Raised if no atoms where added to the lattice. The atoms in the unit cell are needed for computing the neighbors and distances of the lattice.

NoConnectionsError

Raised if no connections have been set up.

Notes

Before calling the analyze function all connections in the lattice have to be set up.

Examples

Construct a square lattice with one atom in the unit cell and nearest neighbors:

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom()
>>> latt.add_connections(num_distances=1)

Call analyze after setting up the connections

>>> latt.analyze()
get_position(nvec=None, alpha=0)[source]

Returns the position for a given translation vector nvec and atom alpha.

Parameters
nvec(N) array_like or int

The translation vector.

alphaint, optional

The atom index, default is 0.

Returns
pos(N) np.ndarray

The position of the transformed lattice site.

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom()
>>> latt.add_connections(1, analyze=True)
>>> latt.get_position([1, 0], alpha=0)
[1. 0.]
get_positions(indices)[source]

Returns the positions for multiple lattice indices.

Parameters
indices(N, D+1) array_like or int

List of lattice indices in the format \((n_1, ..., n_d, α)\).

Returns
pos(N, D) np.ndarray

The positions of the lattice sites.

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom()
>>> latt.add_connections(1)
>>> ind = [[0, 0, 0], [1, 0, 0], [1, 1, 0]]
>>> latt.get_positions(ind)
[[0. 0.]
 [1. 0.]
 [1. 1.]]
estimate_index(pos)[source]

Returns the nearest matching lattice index (n, alpha) for a position.

Parameters
posarray_like or float

The position of the site in world coordinates.

Returns
nvecnp.ndarray

The estimated translation vector \(n\).

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom()
>>> latt.add_connections(1)
>>> latt.estimate_index([1.2, 0.2])
[1 0]
get_neighbors(nvec=None, alpha=0, distidx=0)[source]

Returns the neighour indices of a given site by transforming neighbor data.

Parameters
nvec: (D) array_like or int, optional

The translation vector, the default is the origin.

alpha: int, optional

The atom index, default is 0.

distidx: int, optional

The index of distance to the neighbors, default is 0 (nearest neighbors).

Returns
indices: (N, D+1) np.ndarray

The lattice indices of the neighbor sites.

Raises
NotAnalyzedError

Raised if the lattice distances and base-neighbors haven’t been computed.

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom()
>>> latt.add_connections(1)
>>> latt.get_neighbors(nvec=[0, 0], alpha=0, distidx=0)
[[ 1  0  0]
 [ 0 -1  0]
 [-1  0  0]
 [ 0  1  0]]
get_neighbor_positions(nvec=None, alpha=0, distidx=0)[source]

Returns the neighour positions of a given site by transforming neighbor data.

Parameters
nvec: (D) array_like or int, optional

The translation vector, the default is the origin.

alpha: int, optional

The site index, default is 0.

distidx: int, default

The index of distance to the neighbors, default is 0 (nearest neighbors).

Returns
positions: (N, D) np.ndarray

The positions of the neighbor sites.

Raises
NotAnalyzedError

Raised if the lattice distances and base-neighbors haven’t been computed.

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom()
>>> latt.add_connections(1)
>>> latt.get_neighbor_positions(nvec=[0, 0], alpha=0, distidx=0)
[[ 1.  0.]
 [ 0. -1.]
 [-1.  0.]
 [ 0.  1.]]
get_neighbor_vectors(alpha=0, distidx=0, include_zero=False)[source]

Returns the vectors to the neigbor sites of an atom in the unit cell.

Parameters
alphaint, optional

Index of the base atom. The default is the first atom in the unit cell.

distidxint, default

Index of distance to the neighbors, default is 0 (nearest neighbors).

include_zerobool, optional

Flag if zero-vector is included in result. The default is False.

Returns
vectorsnp.ndarray

The vectors from the site of the atom \(alpha\) to the neighbor sites.

Raises
NotAnalyzedError

Raised if the lattice distances and base-neighbors haven’t been computed.

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom()
>>> latt.add_connections(1)
>>> latt.get_neighbor_vectors(alpha=0, distidx=0)
[[ 1.  0.]
 [ 0. -1.]
 [-1.  0.]
 [ 0.  1.]]
fourier_weights(k, alpha=0, distidx=0)[source]

Returns the Fourier-weight for a given vector.

Parameters
karray_like

The wavevector to compute the lattice Fourier-weights.

alphaint, optional

Index of the base atom. The default is the first atom in the unit cell.

distidxint, default

Index of distance to the neighbors, default is 0 (nearest neighbors).

Returns
weightnp.ndarray
get_base_atom_dict(atleast2d=True)[source]

Returns a dictionary containing the positions for eatch of the base atoms.

Parameters
atleast2dbool, optional

If True, one-dimensional coordinates will be cast to 2D vectors.

Returns
atom_posdict

The positions of the atoms as a dictionary.

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom([0, 0], atom="A")
>>> latt.add_atom([0.5, 0], atom="B")
>>> latt.add_atom([0.5, 0.5], atom="B")
>>> latt.get_base_atom_dict()
{
    Atom(A, radius=0.2, 0): [array([0, 0])],
    Atom(B, radius=0.2, 1): [array([0.5, 0. ]), array([0.5, 0.5])]
}
check_points(points, shape, relative=False, pos=None, tol=0.001)[source]

Returns a mask for the points in the given shape.

Parameters
points: (M, N) np.ndarray

The points in cartesian coordinates.

shape: (N) array_like or int or AbstractShape

shape of finite size lattice to build.

relative: bool, optional

If True the shape will be multiplied by the cell size of the model. The default is True.

pos: (N) array_like or int, optional

Optional position of the section to build. If None the origin is used.

tol: float, optional

The tolerance for checking the points. The default is 1e-3.

Returns
mask: (M) np.ndarray

The mask for the points inside the shape.

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> shape = (2, 2)
>>> points = np.array([[0, 0], [2, 2], [3, 2]])
>>> latt.check_points(points, shape)
[ True  True False]
build_translation_vectors(shape, primitive=False, pos=None, check=True, dtype=None, oversample=0.0)[source]

Constructs the translation vectors \(n\) in a given shape.

Parameters
shape: (N) array_like or int

shape of finite size lattice to build.

primitive: bool, optional

If True the shape will be multiplied by the cell size of the model. The default is True.

pos: (N) array_like or int, optional

Optional position of the section to build. If None the origin is used.

check: bool, optional

If True the positions of the translation vectors are checked and filtered. The default is True. This should only be disabled if filtered later.

dtype: int or np.dtype, optional

Optional data-type for storing the lattice indices. By default, the given limits are checked to determine the smallest possible data-type.

oversample: float, optional

Faktor for upscaling limits for initial index grid. This ensures that all positions are included. Only needed if corner points are missing. The default is 0.

Returns
nvecs: (M, N) np.ndarray

The translation-vectors in lattice-coordinates.

Raises
ValueError

Raised if the dimension of the position doesn’t match the dimension of the lattice.

Examples

>>> latt = LatticeStructure(np.eye(2))
>>> latt.build_translation_vectors((2, 2))
[[0 0]
 [0 1]
 [0 2]
 [1 0]
 [1 1]
 [1 2]
 [2 0]
 [2 1]
 [2 2]]
build_indices(shape, primitive=False, pos=None, check=True, callback=None, dtype=None, return_pos=False)[source]

Constructs the lattice indices .math:(n, α) in the given shape.

Parameters
shape: (N) array_like or int or AbstractShape

shape of finite size lattice to build.

primitive: bool, optional

If True the shape will be multiplied by the cell size of the model. The default is True.

pos: (N) array_like or int, optional

Optional position of the section to build. If None the origin is used.

checkbool, optional

If True the positions of the translation vectors are checked and filtered. The default is True. This should only be disabled if filtered later.

callback: callable, optional

Optional callable for filtering sites. The indices and positions are passed as arguments.

dtype: int or str or np.dtype, optional

Optional data-type for storing the lattice indices. By default, the given limits are checked to determine the smallest possible data-type.

return_pos: bool, optional

Flag if positions should be returned with the indices. This can speed up the building process, since the positions have to be computed here anyway. The default is False.

Returns
indices: (M, N+1) np.ndarray

The lattice indices of the sites in the format .math:(n_1, …, n_d, α).

positions: (M, N) np.ndarray

Corresponding positions. Only returned if return_positions=True.

Raises
ValueError

Raised if the dimension of the position doesn’t match the dimension of the lattice.

Examples

Build indices of a linear chain with two atoms in the unit cell:

>>> latt = LatticeStructure(np.eye(2))
>>> latt.add_atom([0.0, 0.0], "A")
>>> latt.add_atom([0.5, 0.5], "B")
>>> indices, positions = latt.build_indices((2, 1), return_pos=True)

The indices contain the translation vector and the atom index

>>> indices
[[0 0 0]
 [0 0 1]
 [0 1 0]
 [1 0 0]
 [1 0 1]
 [1 1 0]
 [2 0 0]
 [2 1 0]]

The positions are the positions of the atoms in the same order of the indices:

>>> positions
[[0.  0. ]
 [0.5 0.5]
 [0.  1. ]
 [1.  0. ]
 [1.5 0.5]
 [1.  1. ]
 [2.  0. ]
 [2.  1. ]]
compute_neighbors(indices, positions, num_jobs=1)[source]

Computes the neighbors for the given points.

Parameters
indices(N, D+1) array_like

The lattice indices of the sites in the format .math:(n_1, …, n_D, α) where N is the number of sites and D the dimension of the lattice.

positions(N, D) array_like

The positions of the sites in cartesian coordinates where N is the number of sites and D the dimension of the lattice.

num_jobsint, optional

Number of jobs to schedule for parallel processing. If -1 is given all processors are used. The default is 1.

Returns
neighbors: (N, M) np.ndarray

The indices of the neighbors in positions. M is the maximum number of neighbors previously computed in the analyze method.

distances: (N, M) np.ndarray

The corresponding distances of the neighbors.

See also

analyze

Used to pre-compute the base neighbors of the unit cell.

Examples

Construct indices of a one dimensional lattice:

>>> latt = LatticeStructure(1)
>>> latt.add_atom()
>>> latt.add_connections()
>>> indices, positions = latt.build_indices(3, return_pos=True)
>>> positions
[[0.]
 [1.]
 [2.]
 [3.]]

Compute the neighbors of the constructed sites

>>> neighbors, distances = latt.compute_neighbors(indices, positions)
>>> neighbors
[[1 4]
 [2 0]
 [3 1]
 [2 4]]
>>> indices
[[ 1. inf]
 [ 1.  1.]
 [ 1.  1.]
 [ 1. inf]]

The neighbor indices and distances of sites with less than the maximum number of neighbors are filled up with an invalid index (here: 4) and np.inf as distance.

copy()[source]

LatticeStructure : Creates a (deep) copy of the lattice instance.

todict()[source]

Creates a dictionary containing the information of the lattice instance.

Returns
ddict

The information defining the current instance.

classmethod fromdict(d)[source]

Creates a new instance from information stored in a dictionary.

Parameters
ddict

The information defining the current instance.

Returns
lattLatticeStructure

The restored lattice instance.

dumps()[source]

Creates a string containing the information of the lattice instance.

Returns
sstr

The information defining the current instance.

dump(file)[source]

Save the data of the LatticeStructure instance.

Parameters
filestr or int or bytes

File name to store the lattice. If None the hash of the lattice is used.

classmethod load(file)[source]

Load data of a saved LatticeStructure instance.

Parameters
filestr or int or bytes

File name to load the lattice.

Returns
lattLatticeStructure

The lattice restored from the file content.

plot_cell(lw=None, alpha=0.5, cell_lw=None, cell_ls='--', margins=0.1, legend=None, grid=False, show_cell=True, show_vecs=True, show_neighbors=True, con_colors=None, adjustable='box', ax=None, show=False)[source]

Plot the unit cell of the lattice.

Parameters
lwfloat, optional

Line width of the neighbor connections.

alphafloat, optional

The alpha value of the neighbor sites.

cell_lwfloat, optional

The line width used for plotting the unit cell outlines.

cell_lsstr, optional

The line style used for plotting the unit cell outlines.

marginsSequence[float] or float, optional

The margins of the plot.

legendbool, optional

Flag if legend is shown.

gridbool, optional

If True, draw a grid in the plot.

show_neighborsbool, optional

If True the neighbors are plotted.

show_vecsbool, optional

If True the first unit-cell is drawn.

show_cellbool, optional

If True the outlines of the unit cell are plotted.

con_colorsSequence[tuple], optional

list of colors to override the defautl connection color. Each element has to be a tuple with the first two elements being the atom indices of the pair and the third element the color, for example [('A', 'A', 'r')].

adjustableNone or {‘box’, ‘datalim’}, optional

If not None, this defines which parameter will be adjusted to meet the equal aspect ratio. If ‘box’, change the physical dimensions of the Axes. If ‘datalim’, change the x or y data limits. Only applied to 2D plots.

axplt.Axes or plt.Axes3D or None, optional

Parent plot. If None, a new plot is initialized.

showbool, optional

If True, show the resulting plot.