sisl.Grid

class sisl.Grid(shape, bc=None, sc=None, dtype=None, geometry=None)

Bases: SuperCellChild

Real-space grid information with associated geometry.

This grid object handles cell vectors and divisions of said grid.

Parameters
  • shape (float or (3,) of int) – the shape of the grid. A float specifies the grid spacing in Angstrom, while a list of integers specifies the exact grid size.

  • bc (list of int (3, 2) or (3, ), optional) – the boundary conditions for each of the cell’s planes. Default to periodic BC.

  • sc (SuperCell, optional) – the supercell that this grid represents. sc has precedence if both geometry and sc has been specified. Default to [1, 1, 1].

  • dtype (numpy.dtype, optional) – the data-type of the grid, default to numpy.float64.

  • geometry (Geometry, optional) – associated geometry with the grid. If sc has not been passed the supercell will be taken from this geometry.

Examples

>>> grid1 = Grid(0.1, sc=10)
>>> grid2 = Grid(0.1, sc=SuperCell(10))
>>> grid3 = Grid(0.1, sc=SuperCell([10] * 3))
>>> grid1 == grid2
True
>>> grid1 == grid3
True
>>> grid = Grid(0.1, sc=10, dtype=np.complex128)
>>> grid == grid1
False

Methods

add_vacuum(vacuum, axis)

Add vacuum along the axis lattice vector

append(other, axis)

Appends other Grid to this grid along axis

apply(function_, *args, **kwargs)

Applies a function to the grid and returns a new grid.

area(ax0, ax1)

Calculate the area spanned by the two axis ax0 and ax1

average(axis[, weights])

Average grid values along direction axis.

copy([dtype])

Copy the object, possibly changing the data-type

cross_section(idx, axis)

Takes a cross-section of the grid along axis axis

fill(val)

Fill the grid with this value

index(coord[, axis])

Find the grid index for a given coordinate (possibly only along a given lattice vector axis)

index2xyz(index)

Real-space coordinates of indices related to the grid

index_fold(index[, unique])

Converts indices from any placement to only exist in the "primary" grid

index_truncate(index)

Remove indices from outside the grid to only retain indices in the "primary" grid

interp(shape[, order, mode])

Interpolate grid values to a new grid of a different shape

is_orthogonal()

Return true if all cell vectors are linearly independent

isosurface(level[, step_size])

Calculates the isosurface for a given value.

mean(axis[, weights])

Average grid values along direction axis.

mgrid(*slices)

Return a list of indices corresponding to the slices

pyamg_boundary_condition(A, b[, bc])

Attach boundary conditions to the pyamg grid-matrix A with default boundary conditions as specified for this Grid

pyamg_fix(A, b, pyamg_indices, value)

Fix values for the stencil to value.

pyamg_index(index)

Calculate pyamg matrix indices from a list of grid indices

pyamg_source(b, pyamg_indices, value)

Fix the source term to value.

read(sile, *args, **kwargs)

Reads grid from the Sile using read_grid

remove(idx, axis)

Removes certain indices from a specified axis.

remove_part(idx, axis, above)

Removes parts of the grid via above/below designations.

sc_index(*args, **kwargs)

Call local SuperCell object sc_index function

set_bc([boundary, a, b, c])

Set the boundary conditions on the grid

set_boundary([boundary, a, b, c])

Set the boundary conditions on the grid

set_boundary_condition([boundary, a, b, c])

Set the boundary conditions on the grid

set_geom(geometry)

deprecated set_geom

set_geometry(geometry)

Sets the Geometry for the grid.

set_grid(shape[, dtype])

Create the internal grid of certain size.

set_nsc(*args, **kwargs)

Set the number of super-cells in the SuperCell object

set_sc(sc)

Overwrites the local supercell

set_supercell(sc)

Overwrites the local supercell

smooth([r, method, mode])

Make a smoother grid by applying a filter.

sub(idx, axis)

Retains certain indices from a specified axis.

sub_part(idx, axis, above)

Retains parts of the grid via above/below designations.

sum(axis)

Sum grid values along axis axis.

swapaxes(a, b)

Swap two axes in the grid (also swaps axes in the supercell)

tile(reps, axis)

Tile grid to create a bigger one

topyamg([dtype])

Create a pyamg stencil matrix to be used in pyamg

write(sile, *args, **kwargs)

Writes grid to the Sile using write_grid

DIRICHLET

NEUMANN

OPEN

PERIODIC

cell

Returns the inherent SuperCell objects cell

dcell

Voxel cell size

dkind

The data-type of the grid (in str)

dtype

Data-type used in grid

dvolume

Volume of the grid voxel elements

geom

deprecated geometry

icell

Returns the inherent SuperCell objects icell

isc_off

Returns the inherent SuperCell objects isc_off

length

Returns the inherent SuperCell objects length

n_s

Returns the inherent SuperCell objects n_s

nsc

Returns the inherent SuperCell objects nsc

origin

Returns the inherent SuperCell objects origin

origo

Returns the inherent SuperCell objects origin

rcell

Returns the inherent SuperCell objects rcell

sc_off

Returns the inherent SuperCell objects sc_off

shape

Grid shape in \(x\), \(y\), \(z\) directions

size

Total number of elements in the grid

volume

Returns the inherent SuperCell objects volume

DIRICHLET = 3
NEUMANN = 2
OPEN = 4
PERIODIC = 1
__init__(shape, bc=None, sc=None, dtype=None, geometry=None)[source]
add_vacuum(vacuum, axis)

Add vacuum along the axis lattice vector

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

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

append(other, axis)[source]

Appends other Grid to this grid along axis

apply(function_, *args, **kwargs)[source]

Applies a function to the grid and returns a new grid.

You can also apply a function that does not return a grid (maybe you want to do some measurement). In that case, you will get the result instead of a Grid.

Parameters
  • function (str or function) – for a string the full module path to the function should be given. The function that will be called should have the grid as the first argument in its interface.

  • kwargs (args and) – arguments that go directly to the function call

area(ax0, ax1)

Calculate the area spanned by the two axis ax0 and ax1

average(axis, weights=None)[source]

Average grid values along direction axis.

Parameters
  • axis (int) – unit-cell direction to average across

  • weights (array_like, optional) – the weights for the individual axis elements, if boolean it corresponds to 0 and 1 for false/true.

See also

numpy.average

for details regarding the weights argument

property cell

Returns the inherent SuperCell objects cell

copy(dtype=None)[source]

Copy the object, possibly changing the data-type

cross_section(idx, axis)[source]

Takes a cross-section of the grid along axis axis

Remark: This API entry might change to handle arbitrary cuts via rotation of the axis

property dcell

Voxel cell size

property dkind

The data-type of the grid (in str)

property dtype

Data-type used in grid

property dvolume

Volume of the grid voxel elements

fill(val)[source]

Fill the grid with this value

Parameters

val (numpy.dtype) – all grid-points will have this value after execution

property geom

deprecated geometry

property icell

Returns the inherent SuperCell objects icell

index(coord, axis=None)[source]

Find the grid index for a given coordinate (possibly only along a given lattice vector axis)

Parameters
  • coord ((:, 3) or float or Shape) – the coordinate of the axis. If a float is passed axis is also required in which case it corresponds to the length along the lattice vector corresponding to axis. If a Shape a list of coordinates that fits the voxel positions are returned (all internal points also).

  • axis (int, optional) – the axis direction of the index, or for all axes if none.

index2xyz(index)[source]

Real-space coordinates of indices related to the grid

Parameters

index (array_like) – indices for grid-positions

Returns

coordinates of the indices with respect to this grid spacing

Return type

numpy.ndarray

index_fold(index, unique=True)[source]

Converts indices from any placement to only exist in the “primary” grid

Examples

>>> grid = Grid([10, 10, 10])
>>> assert np.all(grid.index_fold([-1, -1, -1]) == 9)
Parameters
  • index (array_like) – indices for grid-positions

  • unique (bool, optional) – if true the returned indices are made unique after having folded the index points

Returns

all indices are then within the shape of the grid

Return type

numpy.ndarray

See also

index_truncate

truncate indices by removing indices outside the primary cell

index_truncate(index)[source]

Remove indices from outside the grid to only retain indices in the “primary” grid

Examples

>>> grid = Grid([10, 10, 10])
>>> assert len(grid.index_truncate([-1, -1, -1])) == 0
Parameters

index (array_like) – indices for grid-positions

Returns

all indices are then within the shape of the grid (others have been removed

Return type

numpy.ndarray

See also

index_fold

fold indices into the primary cell

interp(shape, order=1, mode='wrap', **kwargs)[source]

Interpolate grid values to a new grid of a different shape

It uses the scipy.ndimage.zoom, which creates a finer or more spaced grid using spline interpolation.

Parameters
  • shape (int, array_like of len 3) – the new shape of the grid.

  • order (int 0-5, optional) – the order of the spline interpolation. 1 means linear, 2 quadratic, etc…

  • mode ({'wrap', 'mirror', 'constant', 'reflect', 'nearest'}) – determines how to compute the borders of the grid. The default is 'wrap', which accounts for periodic conditions.

  • **kwargs – optional arguments passed to the interpolation algorithm The interpolation routine is scipy.ndimage.zoom

See also

scipy.ndimage.zoom

method used for interpolation

is_orthogonal()

Return true if all cell vectors are linearly independent

property isc_off

Returns the inherent SuperCell objects isc_off

isosurface(level, step_size=1, **kwargs)[source]

Calculates the isosurface for a given value.

It uses skimage.measure.marching_cubes, so you need to have scikit-image installed.

Parameters
  • level (float) – contour value to search for isosurfaces in the grid. If not given or None, the average of the min and max of the grid is used.

  • step_size (int, optional) – step size in voxels. Larger steps yield faster but coarser results. The result will always be topologically correct though.

  • **kwargs – optional arguments passed directly to skimage.measure.marching_cubes for the calculation of isosurfaces.

Returns

  • numpy array of shape (V, 3) – Verts. Spatial coordinates for V unique mesh vertices.

  • numpy array of shape (n_faces, 3) – Faces. Define triangular faces via referencing vertex indices from verts. This algorithm specifically outputs triangles, so each face has exactly three indices.

  • numpy array of shape (V, 3) – Normals. The normal direction at each vertex, as calculated from the data.

  • numpy array of shape (V, 3) – Values. Gives a measure for the maximum value of the data in the local region near each vertex. This can be used by visualization tools to apply a colormap to the mesh.

See also

skimage.measure.marching_cubes

method used to calculate the isosurface.

property length

Returns the inherent SuperCell objects length

mean(axis, weights=None)

Average grid values along direction axis.

Parameters
  • axis (int) – unit-cell direction to average across

  • weights (array_like, optional) – the weights for the individual axis elements, if boolean it corresponds to 0 and 1 for false/true.

See also

numpy.average

for details regarding the weights argument

classmethod mgrid(*slices)[source]

Return a list of indices corresponding to the slices

The returned values are equivalent to numpy.mgrid but they are returned in a (:, 3) array.

Parameters

*slices (slice or list of int or int) – return a linear list of indices that points to the collective slice made by the passed arguments

Returns

linear indices for each of the sliced values, shape (*, 3)

Return type

numpy.ndarray

property n_s

Returns the inherent SuperCell objects n_s

property nsc

Returns the inherent SuperCell objects nsc

property origin

Returns the inherent SuperCell objects origin

property origo

Returns the inherent SuperCell objects origin

plot

Handles all plotting possibilities for a class

pyamg_boundary_condition(A, b, bc=None)[source]

Attach boundary conditions to the pyamg grid-matrix A with default boundary conditions as specified for this Grid

Parameters
  • A (scipy.sparse.csr_matrix) – sparse matrix describing the grid

  • b (numpy.ndarray) – a vector containing RHS of \(A x = b\) for the solution of the grid stencil

  • bc (list of BC, optional) – the specified boundary conditions. Default to the grid’s boundary conditions, else bc must be a list of elements with elements corresponding to Grid.PERIODIC/Grid.NEUMANN

pyamg_fix(A, b, pyamg_indices, value)[source]

Fix values for the stencil to value.

Parameters
  • A (csr_matrix/csc_matrix) – sparse matrix describing the LHS for the linear system of equations

  • b (numpy.ndarray) – a vector containing RHS of \(A x = b\) for the solution of the grid stencil

  • pyamg_indices (list of int) – the linear pyamg matrix indices where the value of the grid is fixed. I.e. the indices should correspond to returned quantities from pyamg_indices.

  • value (float) – the value of the grid to fix the value at

pyamg_index(index)[source]

Calculate pyamg matrix indices from a list of grid indices

Parameters

index ((:, 3) of int) – a list of indices of the grid along each grid axis

Returns

linear indices for the matrix

Return type

numpy.ndarray

See also

index

query indices from coordinates (directly passable to this method)

mgrid

Grid equivalent to numpy.mgrid. Grid.mgrid returns indices in shapes (:, 3), contrary to numpy’s numpy.mgrid

Raises

ValueError – if any of the passed indices are below 0 or above the number of elements per axis

classmethod pyamg_source(b, pyamg_indices, value)[source]

Fix the source term to value.

Parameters
  • b (numpy.ndarray) – a vector containing RHS of \(A x = b\) for the solution of the grid stencil

  • pyamg_indices (list of int) – the linear pyamg matrix indices where the value of the grid is fixed. I.e. the indices should correspond to returned quantities from pyamg_indices.

property rcell

Returns the inherent SuperCell objects rcell

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

Reads grid from the Sile using read_grid

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

  • * (args passed directly to read_grid(,**)) –

remove(idx, axis)[source]

Removes certain indices from a specified axis.

Works exactly opposite to sub.

Parameters
  • idx (array_like) – the indices of the grid axis axis to be removed

  • axis (int) – the axis segment from which we remove all indices idx

remove_part(idx, axis, above)[source]

Removes parts of the grid via above/below designations.

Works exactly opposite to sub_part

Parameters
  • idx (int) – the index of the grid axis axis to be removed for above=True grid[:idx,…] for above=False grid[idx:,…]

  • axis (int) – the axis segment from which we retain the indices idx

  • above (bool) –

    if True will retain the grid:

    grid[:idx,...]

    else it will retain the grid:

    grid[idx:,...]

sc_index(*args, **kwargs)

Call local SuperCell object sc_index function

property sc_off

Returns the inherent SuperCell objects sc_off

set_bc(boundary=None, a=None, b=None, c=None)[source]

Set the boundary conditions on the grid

Parameters
  • boundary ((3, 2) or (3, ) or int, optional) – boundary condition for all boundaries (or the same for all)

  • a (int or list of int, optional) – boundary condition for the first unit-cell vector direction

  • b (int or list of int, optional) – boundary condition for the second unit-cell vector direction

  • c (int or list of int, optional) – boundary condition for the third unit-cell vector direction

Raises

ValueError – if specifying periodic one one boundary, so must the opposite side.

set_boundary(boundary=None, a=None, b=None, c=None)

Set the boundary conditions on the grid

Parameters
  • boundary ((3, 2) or (3, ) or int, optional) – boundary condition for all boundaries (or the same for all)

  • a (int or list of int, optional) – boundary condition for the first unit-cell vector direction

  • b (int or list of int, optional) – boundary condition for the second unit-cell vector direction

  • c (int or list of int, optional) – boundary condition for the third unit-cell vector direction

Raises

ValueError – if specifying periodic one one boundary, so must the opposite side.

set_boundary_condition(boundary=None, a=None, b=None, c=None)

Set the boundary conditions on the grid

Parameters
  • boundary ((3, 2) or (3, ) or int, optional) – boundary condition for all boundaries (or the same for all)

  • a (int or list of int, optional) – boundary condition for the first unit-cell vector direction

  • b (int or list of int, optional) – boundary condition for the second unit-cell vector direction

  • c (int or list of int, optional) – boundary condition for the third unit-cell vector direction

Raises

ValueError – if specifying periodic one one boundary, so must the opposite side.

set_geom(geometry)[source]

deprecated set_geom

set_geometry(geometry)[source]

Sets the Geometry for the grid.

Setting the Geometry for the grid is a possibility to attach atoms to the grid.

It is not a necessary entity, so passing None is a viable option.

set_grid(shape, dtype=None)[source]

Create the internal grid of certain size.

set_nsc(*args, **kwargs)

Set the number of super-cells in the SuperCell object

See set_nsc for allowed parameters.

See also

SuperCell.set_nsc

the underlying called method

set_sc(sc)

Overwrites the local supercell

set_supercell(sc)

Overwrites the local supercell

property shape

Grid shape in \(x\), \(y\), \(z\) directions

property size

Total number of elements in the grid

smooth(r=0.7, method='gaussian', mode='wrap', **kwargs)[source]

Make a smoother grid by applying a filter.

Parameters
  • r (float or array-like of len 3, optional) –

    the radius of the filter in Angstrom for each axis. If the method is "gaussian", this is the standard deviation!

    If a single float is provided, then the same distance will be used for all axes.

  • method ({'gaussian', 'uniform'}) – the type of filter to apply to smoothen the grid.

  • mode ({'wrap', 'mirror', 'constant', 'reflect', 'nearest'}) – determines how to compute the borders of the grid. The default is wrap, which accounts for periodic conditions.

sub(idx, axis)[source]

Retains certain indices from a specified axis.

Works exactly opposite to remove.

Parameters
  • idx (array_like) – the indices of the grid axis axis to be retained

  • axis (int) – the axis segment from which we retain the indices idx

sub_part(idx, axis, above)[source]

Retains parts of the grid via above/below designations.

Works exactly opposite to remove_part

Parameters
  • idx (int) – the index of the grid axis axis to be retained for above=True grid[idx:,…] for above=False grid[:idx,…]

  • axis (int) – the axis segment from which we retain the indices idx

  • above (bool) –

    if True will retain the grid:

    grid[idx:,...]

    else it will retain the grid:

    grid[:idx,...]

sum(axis)[source]

Sum grid values along axis axis.

Parameters

axis (int) – unit-cell direction to sum across

swapaxes(a, b)[source]

Swap two axes in the grid (also swaps axes in the supercell)

If swapaxes(0,1) it returns the 0 in the 1 values.

Parameters
  • a (int) – axes indices to be swapped

  • b (int) – axes indices to be swapped

tile(reps, axis)[source]

Tile grid to create a bigger one

The atomic indices for the base Geometry will be retained.

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

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

See also

Geometry.tile

equivalent method for Geometry class

topyamg(dtype=None)[source]

Create a pyamg stencil matrix to be used in pyamg

This allows retrieving the grid matrix equivalent of the real-space grid. Subsequently the returned matrix may be used in pyamg for solutions etc.

The pyamg suite is it-self a rather complicated code with many options. For details we refer to pyamg.

Parameters

dtype (numpy.dtype, optional) – data-type used for the sparse matrix, default to use the grid data-type

Returns

  • scipy.sparse.csr_matrix – the stencil for the pyamg solver

  • numpy.ndarray – RHS of the linear system of equations

Examples

This example proves the best method for a variety of cases in regards of the 3D Poisson problem:

>>> grid = Grid(0.01)
>>> A, b = grid.topyamg() # automatically setups the current boundary conditions
>>> # add terms etc. to A and/or b
>>> import pyamg
>>> from scipy.sparse.linalg import cg
>>> ml = pyamg.aggregation.smoothed_aggregation_solver(A, max_levels=1000)
>>> M = ml.aspreconditioner(cycle='W') # pre-conditioner
>>> x, info = cg(A, b, tol=1e-12, M=M)

See also

pyamg_index

convert grid indices into the sparse matrix indices for A

pyamg_fix

fixes stencil for indices and fixes the source for the RHS matrix (uses pyamg_source)

pyamg_source

fix the RHS matrix b to a constant value

pyamg_boundary_condition

setup the sparse matrix A to given boundary conditions (called in this routine)

property volume

Returns the inherent SuperCell objects volume

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

Writes grid to the Sile using write_grid

Parameters

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