Grid

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

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

Attributes

DIRICHLET Constant for defining a Dirichlet boundary condition
NEUMANN Constant for defining a Neumann boundary condition
OPEN Constant for defining an open boundary condition
PERIODIC Constant for defining a periodic boundary condition
cell Returns the inherent SuperCell objects cell
dcell Returns the delta-cell
dkind The data-type of the grid (in str)
dtype Returns the data-type of the grid
dvolume Volume of the grids voxel elements
geom
icell Returns the inherent SuperCell objects icell
isc_off Returns the inherent SuperCell objects isc_off
n_s Returns the inherent SuperCell objects n_s
nsc Returns the inherent SuperCell objects nsc
origo Returns the inherent SuperCell objects origo
rcell Returns the inherent SuperCell objects rcell
sc_off Returns the inherent SuperCell objects sc_off
shape Returns the shape of the grid
size Returns size of the grid
volume Returns the inherent SuperCell objects vol

Methods

__init__(shape[, bc, sc, dtype, geometry]) Initialize self.
add_vacuum(vacuum, axis) Add vacuum along the axis lattice vector
append(other, axis) Appends other Grid to this grid along axis
average(axis[, weights]) Returns the average grid along direction axis.
copy() Returns a copy of the object.
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]) Returns the index along axis axis where coord exists
index2xyz(index) Real-space coordinates of indices related to the grid
interp(shape[, method]) Returns an interpolated version of the grid
is_orthogonal() Return true if all cell vectors are linearly independent
mean(axis[, weights]) Returns the average grid 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) Sets the Geometry for the grid.
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
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) Returns the grid summed along axis axis.
swapaxes(a, b) Returns Grid with swapped axis
topyamg() Create a pyamg stencil matrix to be used in pyamg
write(sile, *args, **kwargs) Writes grid to the Sile using write_grid
DIRICHLET = 3

Constant for defining a Dirichlet boundary condition

NEUMANN = 2

Constant for defining a Neumann boundary condition

OPEN = 4

Constant for defining an open boundary condition

PERIODIC = 1

Constant for defining a periodic boundary condition

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

average(axis, weights=None)[source]

Returns the average grid along direction axis.

Parameters:
axis : int

unit-cell direction to average across

weights : array_like

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
cell

Returns the inherent SuperCell objects cell

copy()[source]

Returns a copy of the object.

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

dcell

Returns the delta-cell

dkind

The data-type of the grid (in str)

dtype

Returns the data-type of the grid

dvolume

Volume of the grids voxel elements

fill(val)[source]

Fill the grid with this value

Parameters:
val : numpy.dtype

all grid-points will have this value after execution

geom
icell

Returns the inherent SuperCell objects icell

index(coord, axis=None)[source]

Returns the index along axis axis where coord exists

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

the axis direction of the index

index2xyz(index)[source]

Real-space coordinates of indices related to the grid

Parameters:
index : array_like

indices for grid-positions

Returns:
numpy.ndarray:

coordinates of the indices with respect to this grid spacing

interp(shape, method='linear', **kwargs)[source]

Returns an interpolated version of the grid

Parameters:
shape : int, array_like

the new shape of the grid

method : str

the method used to perform the interpolation, see scipy.interpolate.interpn for further details.

**kwargs :

optional arguments passed to the interpolation algorithm The interpolation routine is scipy.interpolate.interpn

is_orthogonal()

Return true if all cell vectors are linearly independent

isc_off

Returns the inherent SuperCell objects isc_off

mean(axis, weights=None)

Returns the average grid along direction axis.

Parameters:
axis : int

unit-cell direction to average across

weights : array_like

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:
indices : (:, 3), linear indices for each of the sliced values
n_s

Returns the inherent SuperCell objects n_s

nsc

Returns the inherent SuperCell objects nsc

origo

Returns the inherent SuperCell objects origo

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 : scipy.sparse.csr_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:
pyamg linear indices for the matrix
Raises:
ValueError : if any of the passed indices are below 0 or above the number of elements per axis

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

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

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

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)

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

shape

Returns the shape of the grid

size

Returns size of the grid

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

Raises:
ValueError : if the length of the indices is 0
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]

Returns the grid summed along axis axis.

Parameters:
axis : int

unit-cell direction to sum across

swapaxes(a, b)[source]

Returns Grid with swapped axis

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

topyamg()[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.

Returns:
A : scipy.sparse.csr_matrix which contains the grid stencil for a pyamg solver.
b : numpy.ndarray containing RHS of the linear system of equations.

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)

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

Returns the inherent SuperCell objects vol

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

Writes grid to the Sile using write_grid

Parameters:
sile : Sile, str

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