Grid

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

Object to retain grid information

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

A grid can be periodic and non-periodic.

Attributes

DIRICHLET
NEUMANN
PERIODIC
cell Returns the inherent SuperCell objects cell
dcell Returns the delta-cell
dtype Returns the data-type of the grid
dvol Volume of the grids voxel elements
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
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
vol Returns the inherent SuperCell objects vol

Methods

__init__(shape[, bc, sc, dtype, geom]) Initialize a Grid object.
add_vacuum(vacuum, axis) Add vacuum along the axis lattice vector
append(other, axis) Appends other Grid to this grid along axis
average(axis) 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
index(coord[, axis]) Returns the index along axis axis where coord exists
interp(shape[, method]) Returns an interpolated version of the grid
is_orthogonal() Return true if all cell vectors are linearly independent
mean(axis) Returns the average grid along direction axis
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(geom) Sets the Geometry for the grid.
set_geometry(geom) 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
write(sile, *args, **kwargs) Writes grid to the Sile using write_grid
DIRICHLET = 3
NEUMANN = 2
PERIODIC = 1
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)[source]

Returns the average grid along direction axis

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

dtype

Returns the data-type of the grid

dvol

Volume of the grids voxel elements

index(coord, axis=None)[source]

Returns the index along axis axis where coord exists

Parameters:

coord : array_like or float

the coordinate of the axis

axis : int

the axis direction of the index

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)

Returns the average grid along direction axis

n_s

Returns the inherent SuperCell objects n_s

nsc

Returns the inherent SuperCell objects nsc

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

the indices 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, ) or int, optional

boundary condition for all boundaries (or the same for all)

a: int, optional

boundary condition for the first unit-cell vector direction

b: int, optional

boundary condition for the second unit-cell vector direction

c: int, optional

boundary condition for the third unit-cell vector direction

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

Set the boundary conditions on the grid

Parameters:

boundary: (3, ) or int, optional

boundary condition for all boundaries (or the same for all)

a: int, optional

boundary condition for the first unit-cell vector direction

b: int, optional

boundary condition for the second unit-cell vector direction

c: int, optional

boundary condition for the third unit-cell vector direction

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

Set the boundary conditions on the grid

Parameters:

boundary: (3, ) or int, optional

boundary condition for all boundaries (or the same for all)

a: int, optional

boundary condition for the first unit-cell vector direction

b: int, optional

boundary condition for the second unit-cell vector direction

c: int, optional

boundary condition for the third unit-cell vector direction

set_geom(geom)

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.

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

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

sub_part(idx, axis, above)[source]

Retains parts of the grid via above/below designations.

Works exactly opposite to remove_part

Parameters:

idx : array_like

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

swapaxes(a, b)[source]

Returns Grid with swapped axis

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

vol

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