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

bclist of int (3, 2) or (3, ), optional

the boundary conditions for each of the cell’s planes. Default to periodic BC.

scSuperCell, optional

the supercell that this grid represents. sc has precedence if both geometry and sc has been specified. Default to [1, 1, 1].

dtypenumpy.dtype, optional

the data-type of the grid, default to numpy.float64.

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

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

Short-hand for self.geometry, may be deprecated later

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

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

size

Total number of elements in the grid

volume

Returns the inherent SuperCell objects vol

Methods

__init__(self, shape[, bc, sc, dtype, geometry])

Initialize self.

add_vacuum(self, vacuum, axis)

Add vacuum along the axis lattice vector

append(self, other, axis)

Appends other Grid to this grid along axis

area(self, ax0, ax1)

Calculate the area spanned by the two axis ax0 and ax1

average(self, axis[, weights])

Average grid values along direction axis.

copy(self)

Copy the object

cross_section(self, idx, axis)

Takes a cross-section of the grid along axis axis

fill(self, val)

Fill the grid with this value

index(self, coord[, axis])

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

index2xyz(self, index)

Real-space coordinates of indices related to the grid

index_fold(self, index[, unique])

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

interp(self, shape[, method])

Interpolate grid values to a new grid

is_orthogonal(self)

Return true if all cell vectors are linearly independent

mean(self, axis[, weights])

Average grid values along direction axis.

mgrid(\*slices)

Return a list of indices corresponding to the slices

pyamg_boundary_condition(self, A, b[, bc])

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

pyamg_fix(self, A, b, pyamg_indices, value)

Fix values for the stencil to value.

pyamg_index(self, 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(self, idx, axis)

Removes certain indices from a specified axis.

remove_part(self, idx, axis, above)

Removes parts of the grid via above/below designations.

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

Call local SuperCell object sc_index function

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

Set the boundary conditions on the grid

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

Set the boundary conditions on the grid

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

Set the boundary conditions on the grid

set_geom(self, geometry)

Sets the Geometry for the grid.

set_geometry(self, geometry)

Sets the Geometry for the grid.

set_grid(self, shape[, dtype])

Create the internal grid of certain size.

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

Set the number of super-cells in the SuperCell object

set_sc(self, sc)

Overwrites the local supercell

set_supercell(self, sc)

Overwrites the local supercell

sub(self, idx, axis)

Retains certain indices from a specified axis.

sub_part(self, idx, axis, above)

Retains parts of the grid via above/below designations.

sum(self, axis)

Sum grid values along axis axis.

swapaxes(self, a, b)

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

topyamg(self[, dtype])

Create a pyamg stencil matrix to be used in pyamg

write(self, 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(self, vacuum, axis)

Add vacuum along the axis lattice vector

Parameters
vacuumfloat

amount of vacuum added, in Ang

axisint

the lattice vector to add vacuum along

append(self, other, axis)[source]

Appends other Grid to this grid along axis

area(self, ax0, ax1)

Calculate the area spanned by the two axis ax0 and ax1

average(self, axis, weights=None)[source]

Average grid values along direction axis.

Parameters
axisint

unit-cell direction to average across

weightsarray_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

property cell

Returns the inherent SuperCell objects cell

copy(self)[source]

Copy the object

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

Fill the grid with this value

Parameters
valnumpy.dtype

all grid-points will have this value after execution

property geom

Short-hand for self.geometry, may be deprecated later

property icell

Returns the inherent SuperCell objects icell

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

axisint

the axis direction of the index

index2xyz(self, index)[source]

Real-space coordinates of indices related to the grid

Parameters
indexarray_like

indices for grid-positions

Returns
numpy.ndarray

coordinates of the indices with respect to this grid spacing

index_fold(self, index, unique=True)[source]

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

Parameters
indexarray_like

indices for grid-positions

uniquebool, optional

if true the returned indices are made unique after having folded the index points

Returns
numpy.ndarray

all indices are then within the shape of the grid

Examples

>>> grid = Grid([10, 10, 10])
>>> assert np.all(grid.index2primary([-1, -1, -1]) == 9)
interp(self, shape, method='linear', **kwargs)[source]

Interpolate grid values to a new grid

Parameters
shapeint, array_like

the new shape of the grid

methodstr

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

Return true if all cell vectors are linearly independent

property isc_off

Returns the inherent SuperCell objects isc_off

mean(self, axis, weights=None)

Average grid values along direction axis.

Parameters
axisint

unit-cell direction to average across

weightsarray_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
*slicesslice or list of int or int

return a linear list of indices that points to the collective slice made by the passed arguments

Returns
numpy.ndarray

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

property n_s

Returns the inherent SuperCell objects n_s

property nsc

Returns the inherent SuperCell objects nsc

property origo

Returns the inherent SuperCell objects origo

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

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

Parameters
Ascipy.sparse.csr_matrix

sparse matrix describing the grid

bnumpy.ndarray

a vector containing RHS of \(A x = b\) for the solution of the grid stencil

bclist 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(self, A, b, pyamg_indices, value)[source]

Fix values for the stencil to value.

Parameters
Acsr_matrix/csc_matrix

sparse matrix describing the LHS for the linear system of equations

bnumpy.ndarray

a vector containing RHS of \(A x = b\) for the solution of the grid stencil

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

valuefloat

the value of the grid to fix the value at

pyamg_index(self, 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
numpy.ndarray

linear indices for the matrix

Raises
ValueErrorif 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
bnumpy.ndarray

a vector containing RHS of \(A x = b\) for the solution of the grid stencil

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

Removes certain indices from a specified axis.

Works exactly opposite to sub.

Parameters
idxarray_like

the indices of the grid axis axis to be removed

axisint

the axis segment from which we remove all indices idx

remove_part(self, idx, axis, above)[source]

Removes parts of the grid via above/below designations.

Works exactly opposite to sub_part

Parameters
idxint

the index of the grid axis axis to be removed for above=True grid[:idx,…] for above=False grid[idx:,…]

axisint

the axis segment from which we retain the indices idx

abovebool
if True will retain the grid:

grid[:idx,...]

else it will retain the grid:

grid[idx:,...]

sc_index(self, *args, **kwargs)

Call local SuperCell object sc_index function

property sc_off

Returns the inherent SuperCell objects sc_off

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

aint or list of int, optional

boundary condition for the first unit-cell vector direction

bint or list of int, optional

boundary condition for the second unit-cell vector direction

cint or list of int, optional

boundary condition for the third unit-cell vector direction

Raises
ValueErrorif specifying periodic one one boundary, so must the opposite side.
set_boundary(self, 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)

aint or list of int, optional

boundary condition for the first unit-cell vector direction

bint or list of int, optional

boundary condition for the second unit-cell vector direction

cint or list of int, optional

boundary condition for the third unit-cell vector direction

Raises
ValueErrorif specifying periodic one one boundary, so must the opposite side.
set_boundary_condition(self, 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)

aint or list of int, optional

boundary condition for the first unit-cell vector direction

bint or list of int, optional

boundary condition for the second unit-cell vector direction

cint or list of int, optional

boundary condition for the third unit-cell vector direction

Raises
ValueErrorif specifying periodic one one boundary, so must the opposite side.
set_geom(self, 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(self, 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(self, shape, dtype=None)[source]

Create the internal grid of certain size.

set_nsc(self, *args, **kwargs)

Set the number of super-cells in the SuperCell object

See set_nsc for allowed parameters.

See also

SuperCell.set_nsc

the underlying called method

set_sc(self, sc)

Overwrites the local supercell

set_supercell(self, sc)

Overwrites the local supercell

property shape

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

property size

Total number of elements in the grid

sub(self, idx, axis)[source]

Retains certain indices from a specified axis.

Works exactly opposite to remove.

Parameters
idxarray_like

the indices of the grid axis axis to be retained

axisint

the axis segment from which we retain the indices idx

Raises
ValueErrorif the length of the indices is 0
sub_part(self, idx, axis, above)[source]

Retains parts of the grid via above/below designations.

Works exactly opposite to remove_part

Parameters
idxint

the index of the grid axis axis to be retained for above=True grid[idx:,…] for above=False grid[:idx,…]

axisint

the axis segment from which we retain the indices idx

abovebool
if True will retain the grid:

grid[idx:,...]

else it will retain the grid:

grid[:idx,...]

sum(self, axis)[source]

Sum grid values along axis axis.

Parameters
axisint

unit-cell direction to sum across

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

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

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

Returns the inherent SuperCell objects vol

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

Writes grid to the Sile using write_grid

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