Skip to content

Lattice

Lattice is subclass of NumPy ndarrays that is adapted to represent a numerical field within a discrete 3dimensional space. It adds spatial properties and functionalities to ndimensional arrays such as bounds, unit, neighbourhood assessment and more.

centroids property readonly

Extracts the centroid of cells that have a positive or True value and returns them as a point cloud

Parameters:

Name Type Description Default
threshold float

exclusive lower bound of the values to assume existence of the voxel. Defaults to 0.0.

required

Returns:

Type Description
topogenesis.Cloud

a point cloud representing the centroids of the lattice cells

indices property readonly

Creates one-dimensional integer indices for cells in the lattice

Returns:

Type Description
topogenesis.Lattice

integer lattice of indices

maxbound property readonly

Real maximum bound of the lattice

Returns:

Type Description
numpy.ndarray

real maximum bound

minbound property readonly

Real minimum bound of the lattice

Returns:

Type Description
numpy.ndarray

real minimum bound

apply_stencil(self, stencil, border_condition='pad_outside', padding_value=0)

This method applies the function of a given stencil on the lattice and returns the result in a new lattice.

Parameters:

Name Type Description Default
stencil topogenesis.Stencil

the stencil to be applied on the lattice

required
border_condition str

specifies how the border condition should be treated. The options are {"pad_outside", "pad_inside", "roll"}. "pad_outside" will offset the lattice in every direction by one step, and fill the new cells with the given padding_value and procedes to performing the computation; the resultant lattice in this case has the same shape as the initial lattice. "pad_inside" will perform the computation on the lattice, offsets inside by one cell from each side and returns the remainder cells; the resultant lattice is 2 cell smaller in each dimension than the original lattice. "roll" will assume that the end of each dimension is connected to the beginning of it and interprets the connectivity of the lattice with a rolling approach; the resultant lattice has the same shape is the original lattice. defaults to "pad_outside".

'pad_outside'
padding_value int

value used for padding in case the border_condition is set to "pad_outside".

0

Exceptions:

Type Description
NotImplementedError

"pad_inside" is not implemented yet

Returns:

Type Description
topogenesis.Lattice

a new lattice containing the result of the application of the stencil

Source code in topogenesis/datastructures/datastructures.py
def apply_stencil(self, stencil, border_condition: str = "pad_outside", padding_value: int = 0):
    """This method applies the function of a given stencil on the lattice and returns the result in a new lattice.

    Args:
        stencil (topogenesis.Stencil): 
            the stencil to be applied on the lattice
        border_condition (str, optional): 
            specifies how the border condition should be treated. The options are {"pad_outside", "pad_inside", "roll"}. "pad_outside" will offset the lattice in every direction by one step, and fill the new cells with the given `padding_value` and procedes to performing the computation; the resultant lattice in this case has the same shape as the initial lattice. "pad_inside" will perform the computation on the lattice, offsets inside by one cell from each side and returns the remainder cells; the resultant lattice is 2 cell smaller in each dimension than the original lattice. "roll" will assume that the end of each dimension is connected to the beginning of it and interprets the connectivity of the lattice with a rolling approach; the resultant lattice has the same shape is the original lattice. defaults to "pad_outside".
        padding_value (int, optional): 
            value used for padding in case the `border_condition` is set to "pad_outside".

    Raises:
        NotImplementedError: "pad_inside" is not implemented yet

    Returns:
        topogenesis.Lattice: a new lattice containing the result of the application of the stencil
    """

    if border_condition == "pad_outside":
        # pad the volume with zero in every direction
        padded_arr = np.pad(self, (1, 1),
                            mode='constant',
                            constant_values=(padding_value, padding_value))
        # convert to lattice
        self_padded = to_lattice(padded_arr,
                                 self.minbound - self.unit,
                                 unit=self.unit)

    elif border_condition == "pad_inside":
        raise NotImplementedError

    elif border_condition == "roll":
        self_padded = to_lattice(np.copy(self), self)

    # find the neighbours based on the stencil
    cell_neighbors = self_padded.find_neighbours(stencil)

    # replace neighbours by their value in volume
    self_padded_flat = self_padded.ravel()
    neighbor_values = self_padded_flat[cell_neighbors]

    # apply the function to the neighbour values
    applied = stencil.function(neighbor_values, axis=1)

    # reshape the neighbour applied into the original lattice shape
    applied_3d = applied.reshape(self_padded.shape)

    # reverse the padding procedure
    if border_condition == "pad_outside":
        # trim the padded dimensions
        applied_3d_trimmed = applied_3d[1:-1, 1:-1, 1:-1]

    elif border_condition == "pad_inside":
        raise NotImplementedError

    elif border_condition == "roll":
        applied_3d_trimmed = applied_3d

    return to_lattice(applied_3d_trimmed, self)

arg_apply_stencil(self, arg_lattice, stencil, border_condition='pad_outside', padding_value=0)

Applies the function (should be argument function) of the stencil on the original lattice, extracts the value of the same cell of the argument lattice, fills a new lattice and returns it. If the argument lattice contains one-dimensional ordering of the original lattice, this would function as an argument function.

Parameters:

Name Type Description Default
arg_lattice topogenesis.Lattice

the argument lattice. The values in this lattice will be extracted by the argument function and used to fill the new lattice

required
stencil topogenesis.Stencil

the stencil to be applied on the lattice. This stencil should contain an "argument function"

required
border_condition str

specifies how the border condition should be treated. The options are {"pad_outside", "pad_inside", "roll"}. "pad_outside" will offset the lattice in every direction by one step, and fill the new cells with the given padding_value and procedes to performing the computation; the resultant lattice in this case has the same shape as the initial lattice. "pad_inside" will perform the computation on the lattice, offsets inside by one cell from each side and returns the remainder cells; the resultant lattice is 2 cell smaller in each dimension than the original lattice. "roll" will assume that the end of each dimension is connected to the beginning of it and interprets the connectivity of the lattice with a rolling approach; the resultant lattice has the same shape is the original lattice. defaults to "pad_outside"

'pad_outside'
padding_value int

value used for padding in case the border_condition is set to "pad_outside"

0

Exceptions:

Type Description
ValueError

Main lattice and argument lattice shape should match

NotImplementedError

"pad_inside" is not implemented yet

Returns:

Type Description
topogenesis.Lattice

a new lattice containing the result of the application of the stencil

Source code in topogenesis/datastructures/datastructures.py
def arg_apply_stencil(self, arg_lattice, stencil, border_condition: str = "pad_outside", padding_value: int = 0):
    """Applies the function (should be argument function) of the stencil on the original lattice, extracts the value of the same cell of the argument lattice, fills a new lattice and returns it. If the argument lattice contains one-dimensional ordering of the original lattice, this would function as an argument function. 

    Args:
        arg_lattice (topogenesis.Lattice): 
            the argument lattice. The values in this lattice will be extracted by the argument function and used to fill the new lattice
        stencil (topogenesis.Stencil): 
            the stencil to be applied on the lattice. This stencil should contain an "argument function"
        border_condition (str, optional): 
            specifies how the border condition should be treated. The options are {"pad_outside", "pad_inside", "roll"}. "pad_outside" will offset the lattice in every direction by one step, and fill the new cells with the given `padding_value` and procedes to performing the computation; the resultant lattice in this case has the same shape as the initial lattice. "pad_inside" will perform the computation on the lattice, offsets inside by one cell from each side and returns the remainder cells; the resultant lattice is 2 cell smaller in each dimension than the original lattice. "roll" will assume that the end of each dimension is connected to the beginning of it and interprets the connectivity of the lattice with a rolling approach; the resultant lattice has the same shape is the original lattice. defaults to "pad_outside"
        padding_value (int, optional): 
            value used for padding in case the `border_condition` is set to "pad_outside"

    Raises:
        ValueError: 
            Main lattice and argument lattice shape should match
        NotImplementedError: 
            "pad_inside" is not implemented yet

    Returns:
        topogenesis.Lattice: 
            a new lattice containing the result of the application of the stencil
    """
    if self.shape != arg_lattice.shape:
        raise ValueError(
            "Main lattice and argument lattice shape does not match")

    if border_condition == "pad_outside":
        # pad the volume with padding value in every direction
        padded_arr = np.pad(self, (1, 1),
                            mode='constant',
                            constant_values=(padding_value, padding_value))
        # convert to lattice
        self_padded = to_lattice(padded_arr,
                                 self.minbound - self.unit,
                                 unit=self.unit)
        # pad the argument lattice with padding value in every direction
        padded_arg_arr = np.pad(arg_lattice, (1, 1),
                                mode='constant',
                                constant_values=(padding_value, padding_value))
        # convert to lattice
        arg_lattice_padded = to_lattice(padded_arg_arr,
                                        self.minbound - self.unit,
                                        unit=self.unit)

    elif border_condition == "pad_inside":
        raise NotImplementedError

    elif border_condition == "roll":
        self_padded = to_lattice(np.copy(self), self)
        arg_lattice_padded = to_lattice(np.copy(arg_lattice), arg_lattice)

    # find the neighbours based on the stencil
    cell_neighbors = self_padded.find_neighbours(stencil)

    # replace neighbours by their value in the main lattice
    self_padded_flat = self_padded.ravel()
    neighbor_values = self_padded_flat[cell_neighbors]

    # apply the function to the neighbour values
    applied = stencil.function(neighbor_values, axis=1)
    row_ind = np.arange(applied.size)

    # replace neighbours by their value in the argument latice
    arg_lattice_padded_flat = arg_lattice_padded.ravel()
    arg_neighbor_values = arg_lattice_padded_flat[cell_neighbors]

    # retrieve the values from the argument lattice
    arg_applied = arg_neighbor_values[row_ind, applied]

    # reshape the neighbour applied into the original lattice shape
    arg_applied_3d = arg_applied.reshape(self_padded.shape)

    # reverse the padding procedure
    if border_condition == "pad_outside":
        # trim the padded dimensions
        arg_applied_3d_trimmed = arg_applied_3d[1:-1, 1:-1, 1:-1]

    elif border_condition == "pad_inside":
        raise NotImplementedError

    elif border_condition == "roll":
        arg_applied_3d_trimmed = arg_applied_3d

    return to_lattice(arg_applied_3d_trimmed, self)

boolean_marching_cubes(self)

This is a polygonization method. It converts the lattice to a boolean lattice and runs a boolean marching cube on the lattice.

Returns:

Type Description
topogenesis.Lattice

an integer lattice that contains the tile-id at each cell

Source code in topogenesis/datastructures/datastructures.py
def boolean_marching_cubes(self):
    """This is a polygonization method. It converts the lattice to a boolean lattice and runs a boolean marching cube on the lattice. 

    Returns:
        topogenesis.Lattice: an integer lattice that contains the tile-id at each cell
    """

    # construct the boolean_marching_cubes stencil
    mc_stencil = create_stencil("boolean_marching_cube", 1)

    # retrieve the value of the neighbours
    cell_corners = self.find_neighbours(mc_stencil, order="F")

    # converting volume value (TODO: this needs to become a method of its own)
    volume_flat = self.ravel()
    volume_flat[volume_flat > 0.0] = 1
    volume_flat[volume_flat <= 0.0] = 0

    # replace neighbours by their value in volume
    neighbor_values = volume_flat[cell_corners]

    # computing the cell tile id
    # the powers of 2 in an array
    legend = 2**np.arange(8)

    # multiply the corner with the power of two, sum them, and reshape to the original volume shape
    tile_id = np.sum(legend * neighbor_values,
                     axis=1).reshape(self.shape)

    # drop the last column, row and page (since cube-grid is 1 less than the voxel grid in every dimension)
    # TODO consider that by implementing the origin attribute in lattice this may have to change
    cube_grid = tile_id[:-1, :-1, :-1]

    # convert the array to lattice
    cube_lattice = to_lattice(
        cube_grid, minbound=self.minbound, unit=self.unit)

    return cube_lattice

centroids_threshold(self, threshold=0.0)

Extracts the centroid of cells that have a positive or True value and returns them as a point cloud

Parameters:

Name Type Description Default
threshold float

exclusive lower bound of the values to assume existence of the voxel. Defaults to 0.0.

0.0

Returns:

Type Description
topogenesis.Cloud

a point cloud representing the centroids of the lattice cells

Source code in topogenesis/datastructures/datastructures.py
def centroids_threshold(self, threshold=0.0):
    """Extracts the centroid of cells that have a positive or True value and returns them as a point cloud

    Args:
        threshold (float, optional): exclusive lower bound of the values to assume existence of the voxel. Defaults to 0.0.

    Returns:
        topogenesis.Cloud: a point cloud representing the centroids of the lattice cells
    """

    # extract the indices of the True values # with sparse matrix we don't need to search
    point_array = np.argwhere(self > threshold)
    # convert to float
    point_array = point_array.astype(float)
    # scale by unit
    point_array *= self.unit
    # translate by minimum
    point_array += self.minbound
    # orient the points
    if np.sum(np.abs(self.orient[:3])) > 0.001 or self.orient[3] != 1: # check if the orientation is required
        r = sp_rotation.from_quat(self.orient)
        point_array = r.apply(point_array)
    # return as a point cloud
    return cloud(point_array, dtype=float)

fast_notebook_vis(self, plot, show_outline=True, show_centroids=True)

Adds the basic lattice features to a pyvista ITK plotter and returns it. ITK plotters are specifically used in notebooks to plot the geometry inside the notebook environment It is mainly used to rapidly visualize the content of the lattice for visual confirmation

Parameters:

Name Type Description Default
plot pyvista.PlotterITK

a pyvista ITK plotter

required
show_outline bool

If True, adds the bounding box of the lattice to the plot

True
show_centroids bool

If True, adds the centroid of cells to the plot

True

Returns:

Type Description
pyvista.PlotterITK

pyvista ITK plotter containing lattice features such as cells, bounding box, cell centroids

Usage Example:

p = pyvista.PlotterITK()
lattice.fast_notebook_vis(p)    
Source code in topogenesis/datastructures/datastructures.py
def fast_notebook_vis(self, plot, show_outline: bool = True, show_centroids: bool = True):
    """Adds the basic lattice features to a pyvista ITK plotter and returns it. 
    ITK plotters are specifically used in notebooks to plot the geometry inside 
    the notebook environment It is mainly used to rapidly visualize the content 
    of the lattice for visual confirmation

    Args:
        plot (pyvista.PlotterITK): a pyvista ITK plotter
        show_outline (bool, optional): If `True`, adds the bounding box of the lattice to the plot
        show_centroids (bool, optional): If `True`, adds the centroid of cells to the plot

    Returns:
        pyvista.PlotterITK: pyvista ITK plotter containing lattice features such as cells, bounding box, cell centroids

    ** Usage Example: **
    ```python
    p = pyvista.PlotterITK()
    lattice.fast_notebook_vis(p)    
    ```
    """

    # Set the grid dimensions: shape + 1 because we want to inject our values on the CELL data
    grid = pv.UniformGrid()
    grid.dimensions = np.array(self.shape) + 1
    # The bottom left corner of the data set
    grid.origin = self.minbound - self.unit * 0.5
    grid.spacing = self.unit  # These are the cell sizes along each axis
    # Add the data values to the cell data
    grid.cell_arrays["values"] = self.flatten(
        order="F").astype(float)  # Flatten the array!
    # filtering the voxels
    threshed = grid.threshold([0.9, 1.1])

    # adding the voxels: light red
    plot.add_mesh(threshed, color="#ff8fa3", opacity=0.3)
    # plot.add_mesh(threshed, show_edges=True, color="#ff8fa3", opacity=0.3, label="Cells")

    if show_outline:
        # adding the boundingbox wireframe
        plot.add_mesh(grid.outline(), color="grey")
        # plot.add_mesh(grid.outline(), color="grey", label="Domain")

    if show_centroids:
        # adding the voxel centroids: red
        plot.add_points(pv.PolyData(self.centroids), color='#ff244c')
        # plot.add_mesh(pv.PolyData(self.centroids), color='#ff244c', point_size=5, render_points_as_spheres=True, label="Cell Centroids")

    return plot

fast_vis(self, plot, show_outline=True, show_centroids=True, color='#ff8fa3', opacity=0.3)

Adds the basic lattice features to a pyvista plotter and returns it. It is mainly used to rapidly visualize the content of the lattice for visual confirmation

Parameters:

Name Type Description Default
plot pyvista.Plotter

a pyvista plotter

required
show_outline bool

If True, adds the bounding box of the lattice to the plot

True
show_centroids bool

If True, adds the centroid of cells to the plot

True

Returns:

Type Description
pyvista.Plotter

the same pyvista plotter containing lattice features

Usage Example:

p = pyvista.Plotter()
lattice.fast_vis(p)
Source code in topogenesis/datastructures/datastructures.py
def fast_vis(self, plot, show_outline: bool = True, show_centroids: bool = True, color = "#ff8fa3", opacity=0.3):
    """Adds the basic lattice features to a pyvista plotter and returns it. 
    It is mainly used to rapidly visualize the content of the lattice 
    for visual confirmation

    Args:
        plot (pyvista.Plotter): a pyvista plotter
        show_outline (bool, optional): If `True`, adds the bounding box of the lattice to the plot
        show_centroids (bool, optional): If `True`, adds the centroid of cells to the plot

    Returns:
        pyvista.Plotter: the same pyvista plotter containing lattice features

    ** Usage Example: **
    ```python
    p = pyvista.Plotter()
    lattice.fast_vis(p)
    ```
    """
    # TODO: Add documentation for color and opacity

    # Set the grid dimensions: shape + 1 because we want to inject our values on the CELL data
    grid = pv.UniformGrid()
    grid.dimensions = np.array(self.shape) + 1
    # The bottom left corner of the data set
    grid.origin = self.minbound - self.unit * 0.5
    grid.spacing = self.unit  # These are the cell sizes along each axis
    # Add the data values to the cell data
    grid.cell_arrays["values"] = self.flatten(
        order="F").astype(float)  # Flatten the array!
    # filtering the voxels
    threshed = grid.threshold([0.9, 1.1])

    # applying the orientation of the lattice
    if np.sum(np.abs(self.orient[:3])) > 0.001 or self.orient[3] != 1: # check if the orientation is required
        Rz = sp_rotation.from_quat(self.orient)
        threshed.points = Rz.apply(threshed.points)

    # adding the voxels: light red
    plot.add_mesh(threshed, show_edges=True, color=color,
                  opacity=opacity, label="Cells")

    if show_outline:
        # adding the boundingbox wireframe
        wireframe = grid.outline()
        Rz = sp_rotation.from_quat(self.orient)
        wireframe.points = Rz.apply(wireframe.points)
        plot.add_mesh(wireframe, color="grey", label="Domain")

    if show_centroids:
        # adding the voxel centroids: red
        plot.add_mesh(pv.PolyData(self.centroids), color='#ff244c', point_size=5,
                      render_points_as_spheres=True, label="Cell Centroids")

    return plot

find_neighbours(self, stencil, order='dist')

Given an stencil, this method will return the neighbours of all cells with regard to the stencil specification in a numpy 2D array with each row corresponding to one cell in lattice.

Parameters:

Name Type Description Default
stencil topogenesis.Stencil

Stencil that describes the neighbourhood

required
order str

the order of neighbours is one of {"dist", "C", "F"}. 'dist' sorts the neighbours based on the distance from origin cell, ‘C’ sorts in row-major (C-style) order, and ‘F’ sorts in column-major (Fortran- style) order. defaults to "dist"

'dist'

Returns:

Type Description
numpy.ndarray

2D array describing the neighbours of each cell in the row

Source code in topogenesis/datastructures/datastructures.py
def find_neighbours(self, stencil, order: str = "dist"):
    """Given an stencil, this method will return the neighbours of all cells with regard to the stencil specification in a numpy 2D array with each row corresponding to one cell in lattice.

    Args:
        stencil (topogenesis.Stencil): Stencil that describes the neighbourhood
        order (str, optional): the order of neighbours is one of {"dist", "C", "F"}. 'dist' sorts the neighbours based on the distance from origin cell, ‘C’ sorts in row-major (C-style) order, and ‘F’ sorts in column-major (Fortran- style) order. defaults to "dist"

    Returns:
        numpy.ndarray:  2D array describing the neighbours of each cell in the row
    """
    # TODO: kwarg for returning 1d or 3d indices
    # the id of voxels (0,1,2, ... n)
    self_ind = self.indices

    # calculating all the possible shifts to apply to the array
    shifts = stencil.expand(order)

    # gathering all the replacements in the columns
    replaced_columns = [
        np.roll(self_ind, shift, np.arange(3)).ravel() for shift in shifts]

    # stacking the columns
    cell_neighbors = np.stack(replaced_columns, axis=-1)

    return cell_neighbors

find_neighbours_masked(self, stencil, loc, order='dist', mask=None, border_condition='standard', id_type='1D')

Given an stencil, this method will return the neighbours of one specific cell (specified with loc) with regard to the neighbourhood that stencil defines. The neighbours can be returned with 1D or 3D indices

Parameters:

Name Type Description Default
stencil topogenesis.Stencil

Stencil that describes the neighbourhood

required
loc numpy.ndarray

The location (3D index) of the cells that it's neighbours are desired.

required
order str

the order of neighbours is one of {"dist", "C", "F"}. 'dist' sorts the neighbours based on the distance from origin cell, ‘C’ sorts in row-major (C-style) order, and ‘F’ sorts in column-major (Fortran- style) order. defaults to "dist"

'dist'
mask numpy.ndarray

Not implemented yet. Defaults to None.

None
border_condition str

specifies how the border condition should be treated. The options are {"standard", "roll"}. "standard" will assume that the cells at the border of bound of the lattice and they have less neighbours compared to the cells in the middle of the lattice. "roll" will assume that the end of each dimension is connected to the beginning of it and interprets the connectivity of the lattice with a rolling approach; the resultant lattice has the same shape is the original lattice. defaults to "standard"

'standard'
id_type str

specifies how the neighbours are specified and returened. The options are {"1D", "3D"}. "1D" will return a one-dimensional index of the cell, similar to the index that is assigned to each cell when lattice.indices is called. "3D" will return the three-dimensional index of the cell which is the integer coordinates of the cell as well. Defaults to "1D".

'1D'

Exceptions:

Type Description
NotImplementedError

in case that "roll" option is chosen for "border_condition"

NotImplementedError

in case that "mask" option is used

Returns:

Type Description
numpy.ndarray

array of the indices of the neighbours.

Source code in topogenesis/datastructures/datastructures.py
def find_neighbours_masked(self, stencil, loc, order="dist", mask=None, border_condition="standard", id_type="1D"):
    """Given an stencil, this method will return the neighbours of one specific cell (specified with loc) with regard to the neighbourhood that stencil defines. The neighbours can be returned with 1D or 3D indices

    Args:
        stencil (topogenesis.Stencil): Stencil that describes the neighbourhood
        loc (numpy.ndarray): The location (3D index) of the cells that it's neighbours are desired.
        order (str, optional): the order of neighbours is one of {"dist", "C", "F"}. 'dist' sorts the neighbours based on the distance from origin cell, ‘C’ sorts in row-major (C-style) order, and ‘F’ sorts in column-major (Fortran- style) order. defaults to "dist"
        mask (numpy.ndarray, optional): Not implemented yet. Defaults to None.
        border_condition (str, optional): 
            specifies how the border condition should be treated. The options are {"standard", "roll"}. "standard" will assume that the cells at the border of bound of the lattice and they have less neighbours compared to the cells in the middle of the lattice. "roll" will assume that the end of each dimension is connected to the beginning of it and interprets the connectivity of the lattice with a rolling approach; the resultant lattice has the same shape is the original lattice. defaults to "standard"
        id_type (str, optional): specifies how the neighbours are specified and returened. The options are {"1D", "3D"}. "1D" will return a one-dimensional index of the cell, similar to the index that is assigned to each cell when lattice.indices is called. "3D" will return the three-dimensional index of the cell which is the integer coordinates of the cell as well. Defaults to "1D".

    Raises:
        NotImplementedError: in case that "roll" option is chosen for "border_condition"
        NotImplementedError: in case that "mask" option is used

    Returns:
        numpy.ndarray: array of the indices of the neighbours.
    """
    # the id of voxels (0,1,2, ... n)
    self_ind = self.indices

    # find the bounds of window around the specified location
    win_min = loc + stencil.minbound
    win_max = loc + stencil.maxbound + 1
    if border_condition == "standard":
        # ensure win_min is not less than zero
        win_min = np.maximum(win_min, np.array([0, 0, 0]))
        # ensure win_max is not more than shape -1
        win_max = np.minimum(win_max, np.array(self_ind.shape))

    # TODO
    if border_condition == "roll":
        raise NotImplementedError

    # TODO
    if mask != None:
        raise NotImplementedError

    self_ind = self_ind[win_min[0]: win_max[0],
                        win_min[1]: win_max[1],
                        win_min[2]: win_max[2]]

    # TODO: Sometimes, self_ind turns to be zero. need to be debugged
    # find the new 1D-index of the location
    new_loc_ind = np.ravel_multi_index(
        tuple(stencil.origin), self_ind.shape)

    # calculating all the possible shifts to apply to the array
    shifts = stencil.expand(order)

    # gathering all the replacements in the columns
    replaced_columns = [
        np.roll(self_ind, shift, np.arange(3)).ravel() for shift in shifts]

    # stacking the columns
    cell_neighbors = np.stack(replaced_columns, axis=-1)

    # extract 1D ids
    neighs_1d_id = cell_neighbors[new_loc_ind]

    if id_type == "1D":
        return neighs_1d_id

    elif id_type == "3D":
        # convert 1D index to 3D index
        neigh_3d_id = np.array(
            np.unravel_index(neighs_1d_id, self.shape)).T
        return neigh_3d_id

to_csv(self, filepath)

This method saves the lattice to a csv file

Parameters:

Name Type Description Default
filepath str

path to the csv file

required
Source code in topogenesis/datastructures/datastructures.py
def to_csv(self, filepath: str):
    """This method saves the lattice to a csv file

    Args:
        filepath: path to the csv file
    """
    # volume to panda dataframe
    vol_df = self.to_pandas()

    # specifying metadata and transposing it
    metadata = pd.DataFrame({
        'minbound': self.minbound,
        'shape': np.array(self.shape),
        'unit': self.unit,
    })

    with open(filepath, 'w') as df_out:

        metadata.to_csv(df_out, index=False,
                        header=True, float_format='%g', line_terminator='\n')

        df_out.write('\n')

        vol_df.to_csv(df_out,
                      index=False,
                      float_format='%g', line_terminator='\n')

to_pandas(self)

This methods returns a pandas dataframe containing the lattice information with integer indices and value of the cell as columns.

Returns:

Type Description
pandas.Dataframe

lattice represented in pandas dataframe

Source code in topogenesis/datastructures/datastructures.py
def to_pandas(self):
    """This methods returns a pandas dataframe containing the lattice information with integer indices and value of the cell as columns.

    Returns:
        pandas.Dataframe: lattice represented in pandas dataframe
    """
    # get the indices of the voxels
    vol_3d_ind = np.indices(self.shape)

    # flatten except the last dimension
    vol_3d_ind_flat = vol_3d_ind.transpose(1, 2, 3, 0).reshape(-1, 3)

    # flatten the volume
    vol_flat = self.ravel()

    # volume data to panda dataframe
    vol_df = pd.DataFrame(
        {
            'IX': vol_3d_ind_flat[:, 0],
            'IY': vol_3d_ind_flat[:, 1],
            'IZ': vol_3d_ind_flat[:, 2],
            'value': vol_flat,
        })
    return vol_df

topoGenesis DataStructure

to_lattice(a, minbound, unit=1, orient=array([0., 0., 0., 1.]))

Converts a numpy array into a lattice

Parameters:

Name Type Description Default
a numpy.ndarray

array

required
minbound numpy.ndarray

describing the minimum bound of the lattice in continuous space

required
unit int

the unit size of the lattice

1

Returns:

Type Description
topogenesis.Lattice

the lattice representation of the array

Source code in topogenesis/datastructures/datastructures.py
def to_lattice(a, minbound: np.ndarray, unit=1, orient=np.array([0., 0., 0., 1.])) -> lattice:
    """Converts a numpy array into a lattice

    Args:
        a (numpy.ndarray): 
            array
        minbound (numpy.ndarray): 
            describing the minimum bound of the lattice in continuous space
        unit (int, optional): 
            the unit size of the lattice

    Returns:
        topogenesis.Lattice: the lattice representation of the array
    """
    # check if the minbound is a lattice
    if type(minbound) is lattice:
        l = minbound
        unit = l.unit
        minbound = l.minbound
        orient = l.orient

    # check if minbound is an np array
    elif type(minbound) is not np.ndarray:
        minbound = np.array(minbound)

    # compute the bounds based on the minbound
    bounds = np.array([minbound, minbound + unit * (np.array(a.shape) - 1)])

    # construct a lattice
    array_lattice = lattice(bounds, unit=unit, dtype=a.dtype, orient=orient)
    array_lattice[:, :, :] = a[:, :, :]
    return array_lattice