|  | 
|  | __init__ (self, net=None, dims=None) | 
|  | emptyLatticeMask (self) | 
|  | ellipsoid (self, xs=None, radius=None, angles=None, center=None) | 
|  | cylinder (self, radius, length, xs=None, angles=None, center=None) | 
|  | spoke (self, x0, length, spoke_radius, r, phi, theta) | 
|  | capsule (self, length, width, xs=None, angles=None, center=None) | 
|  | box (self, lx, ly, lz, xs=None, angles=None, center=None) | 
|  | compose (self, *siteSpec, net=None) | 
|  | se6 (self) | 
|  | se26 (self) | 
|  | dilate (cls, binaryMask, radius=None, se=None) | 
|  | erode (cls, binaryMask, radius=None, se=None) | 
|  | closing (cls, binaryMask, radius=None, se=None, radius1=None, se1=None) | 
|  | opening (cls, binaryMask, radius=None, se=None, radius1=None, se1=None) | 
|  | convexHull (cls, binaryMask) | 
|  | 
|  | showBinaryLattices (binLattices, manualColor=None, filterFunctions=None, mode="widget") | 
|  | showStack (binLattices, plane='xz', scl=None, maxWidth=600, maxHeight=600) | 
|  | transformGrid (xs, x0, alpha, beta, gamma) | 
|  | octoStructElem (r) | 
|  | sphereStructElem (r) | 
|  | 
|  | nx | 
|  | ny | 
|  | nz = net.siteLattice.shape | 
|  | shape = net.siteLattice.shape | 
|  | xs = np.mgrid[0:self.nx, 0:self.ny, 0:self.nz] | 
|  | center = np.array( [self.nx//2, self.ny//2, self.nz//2]) | 
|  | origin = np.zeros( 3 ) | 
|  | 
|  | _parseArgs (self, angles, center, xs) | 
|  | _morphApply (cls, bi, radius, se, fn) | 
|  | _morphApplyTwice (cls, bi, radius1, se1, fn1, radius2, se2, fn2) | 
Helper object to design site lattice geometry
 ◆ __init__()
      
        
          | jLM.RegionBuilder.RegionBuilder.__init__ | ( |  | self, | 
        
          |  |  |  | net = None, | 
        
          |  |  |  | dims = None ) | 
      
 
Helper object to design site lattice geometry
Args:
    net (:py:class:`~jLM.RDME.SpatialModel`): 
        If present, take dimensions from simulation 
    dims ([int]): 
        Dimensions of lattice
 
 
 
◆ _morphApply()
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder._morphApply | ( |  | cls, |  
          |  |  |  | bi, |  
          |  |  |  | radius, |  
          |  |  |  | se, |  
          |  |  |  | fn ) |  | protected | 
 
 
◆ _morphApplyTwice()
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder._morphApplyTwice | ( |  | cls, |  
          |  |  |  | bi, |  
          |  |  |  | radius1, |  
          |  |  |  | se1, |  
          |  |  |  | fn1, |  
          |  |  |  | radius2, |  
          |  |  |  | se2, |  
          |  |  |  | fn2 ) |  | protected | 
 
 
◆ _parseArgs()
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder._parseArgs | ( |  | self, |  
          |  |  |  | angles, |  
          |  |  |  | center, |  
          |  |  |  | xs ) |  | protected | 
 
 
◆ box()
      
        
          | jLM.RegionBuilder.RegionBuilder.box | ( |  | self, | 
        
          |  |  |  | lx, | 
        
          |  |  |  | ly, | 
        
          |  |  |  | lz, | 
        
          |  |  |  | xs = None, | 
        
          |  |  |  | angles = None, | 
        
          |  |  |  | center = None ) | 
      
 
Construct a rectangular cuboid mask
Args:
    lx (float):
        Length in x-direction
    ly (float):
        Length in y-direction
    lz (float):
        Length in z-direction
    center (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        Centroid
    angles (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        [alpha, beta, gamma] ZXZ Euler angles
    xs (:py:class:`~numpy.ndarray(shape=(3,nx,ny,nz), dtype=int)`):
        Index grid of dimensions [3, nx, ny, nz]
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`:
        Lattice mask
 
 
 
◆ capsule()
      
        
          | jLM.RegionBuilder.RegionBuilder.capsule | ( |  | self, | 
        
          |  |  |  | length, | 
        
          |  |  |  | width, | 
        
          |  |  |  | xs = None, | 
        
          |  |  |  | angles = None, | 
        
          |  |  |  | center = None ) | 
      
 
Construct spherocylinder mask
Args:
    length (float):
        Length including endcaps
    width (float):
        Diameter of cylindrical region
    center (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        Centroid
    angles (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        [alpha, beta, gamma] ZXZ' Euler angles
    xs (:py:class:`~numpy.ndarray(shape=(3,nx,ny,nz), dtype=int)`): 
        Index grid of dimensions [3, nx, ny, nz]
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Lattice mask
 
 
 
◆ closing()
      
        
          | jLM.RegionBuilder.RegionBuilder.closing | ( |  | cls, | 
        
          |  |  |  | binaryMask, | 
        
          |  |  |  | radius = None, | 
        
          |  |  |  | se = None, | 
        
          |  |  |  | radius1 = None, | 
        
          |  |  |  | se1 = None ) | 
      
 
Morphological closing of a binary mask
Args:
    binaryMask (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Binary mask
    radius (int):
        If provided use a 6-connected structuring element iterated 
        `radius` times.
    radius1 (int):
        If provided use the iterated SE for the erosion (optional)
    se (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Structuring element
    se1 (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Structuring element for erosion (optional)
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Closed lattice
 
 
 
◆ compose()
      
        
          | jLM.RegionBuilder.RegionBuilder.compose | ( |  | self, | 
        
          |  |  | * | siteSpec, | 
        
          |  |  |  | net = None ) | 
      
 
Compose a series of binary masks into a site lattice
This function takes an indefinite number of (region, mask) tuples. The 
lattice is created by setting the lattice to the index of region over 
all masked subvolumes in the order given.  
Args:
    *siteSpec ((:py:class:`~jLM.Types.Region`, :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`)): 
        Region, mask tuple
    net (:py:class:`~jLM.RDME.SpatialModel`): 
        If given, modify this model's site lattice
 
 
 
◆ convexHull()
      
        
          | jLM.RegionBuilder.RegionBuilder.convexHull | ( |  | cls, | 
        
          |  |  |  | binaryMask ) | 
      
 
Convex hull of lattice
Args:
    binaryMask (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Binary mask
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Convex hull of lattice
 
 
 
◆ cylinder()
      
        
          | jLM.RegionBuilder.RegionBuilder.cylinder | ( |  | self, | 
        
          |  |  |  | radius, | 
        
          |  |  |  | length, | 
        
          |  |  |  | xs = None, | 
        
          |  |  |  | angles = None, | 
        
          |  |  |  | center = None ) | 
      
 
Construct cylinder mask
Args:
    radius (float):
        Radius
    length (float):
        Length
    center (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        Centroid
    angles (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        [alpha, beta, gamma] ZXZ' Euler angles
    xs (:py:class:`~numpy.ndarray(shape=(3,nx,ny,nz), dtype=int)`): 
        Index grid of dimensions [3, nx, ny, nz]
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Lattice mask
 
 
 
◆ dilate()
      
        
          | jLM.RegionBuilder.RegionBuilder.dilate | ( |  | cls, | 
        
          |  |  |  | binaryMask, | 
        
          |  |  |  | radius = None, | 
        
          |  |  |  | se = None ) | 
      
 
Morphological dialation of a binary mask
Args:
    binaryMask (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Binary mask
    radius (int):   
        If provided use a 6-connected structuring element iterated `radius` times.
    se (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Structuring element
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Dilated lattice
 
 
 
◆ ellipsoid()
      
        
          | jLM.RegionBuilder.RegionBuilder.ellipsoid | ( |  | self, | 
        
          |  |  |  | xs = None, | 
        
          |  |  |  | radius = None, | 
        
          |  |  |  | angles = None, | 
        
          |  |  |  | center = None ) | 
      
 
Construct ellipsoid mask
Args:
    radius (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`): 
        Semiaxes of ellipse
    center (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        Centroid
    angles (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):   
        [alpha, beta, gamma] ZXZ' Euler angles
    xs (:py:class:`~numpy.ndarray(shape=(3,nx,ny,nz), dtype=int)`): 
        Index grid of dimensions [3, nx, ny, nz]
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Lattice mask
 
 
 
◆ emptyLatticeMask()
      
        
          | jLM.RegionBuilder.RegionBuilder.emptyLatticeMask | ( |  | self | ) |  | 
      
 
Return an empty lattice mask
 
 
 
◆ erode()
      
        
          | jLM.RegionBuilder.RegionBuilder.erode | ( |  | cls, | 
        
          |  |  |  | binaryMask, | 
        
          |  |  |  | radius = None, | 
        
          |  |  |  | se = None ) | 
      
 
Morphological erosion of a binary mask
Args:
    binaryMask (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Binary mask
    radius (int):
        If provided use a 6-connected structuring element iterated 
        `radius` times.
    se (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Structuring element
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Eroded lattice
 
 
 
◆ octoStructElem()
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder.octoStructElem | ( |  | r | ) |  |  | static | 
 
Iterated 6-connected structuring element
Args:
    r (int):
        Number of iterations
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Structuring element
 
 
 
◆ opening()
      
        
          | jLM.RegionBuilder.RegionBuilder.opening | ( |  | cls, | 
        
          |  |  |  | binaryMask, | 
        
          |  |  |  | radius = None, | 
        
          |  |  |  | se = None, | 
        
          |  |  |  | radius1 = None, | 
        
          |  |  |  | se1 = None ) | 
      
 
Morphological opening of a binary mask
Args:
    binaryMask (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Binary mask
    radius (int):   
        If provided use a 6-connected structuring element iterated 
        `radius` times.
    radius1 (int):  
        If provided use the iterated SE for the dilation (optional)
    se (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Structuring element
    se1 (:py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`):
        Structuring element for dilation (optional)
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`:
        Opened lattice
 
 
 
◆ se26()
      
        
          | jLM.RegionBuilder.RegionBuilder.se26 | ( |  | self | ) |  | 
      
 
Structuring element connecting all 26 neighbors
 
 
 
◆ se6()
      
        
          | jLM.RegionBuilder.RegionBuilder.se6 | ( |  | self | ) |  | 
      
 
Structuring element connecting all 6 nearest neighbors
 
 
 
◆ showBinaryLattices()
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder.showBinaryLattices | ( |  | binLattices, |  
          |  |  |  | manualColor = None, |  
          |  |  |  | filterFunctions = None, |  
          |  |  |  | mode = "widget" ) |  | static | 
 
3-D lattice mask viewer
Lattices can be given as a single binary mask, a list of masks, or a 
dict of masks.  The display mode can be "widget", which displays in the 
notebook, "download_x3d", which opens a download link in the notebook to
the X3D scene, or "download_html", which opens a download link in the 
notebook to a standalone HTML file. 
Partial lattices can be shown by applying a filter function.  To hide 
parts of the lattice, `filterFunctions` can be specified.  This option 
takes a list of functions which map from a (x,y,z) mesh grid to 
a [nx,ny,nz] boolean mask where only subvolumes marked True are shown. 
To only show volumes whose z coordinate are less than 32, the function
>>> def zfilter(x,y,z):
>>>     return z<32
is used. Here the arguments x,y,z are of type :py:class:`numpy.ndarray` 
and a boolean lattice is returned.  The filter functions are specified 
using a dictionary where the keys correspond to the names of the 
lattices given. These names come from the dictionary if a dict of 
lattices is given, `latticeXY` if a list of lattices is given (XY
being a zero-padded index into the list), or `lattice` if a single
mask is given.
Args:
    binLattices:
        Binary lattice or collection of binary lattices
    filterFunctions (dict(str=func)):
        Dict of functions indexed by name of lattice
    mode (str):
        View mode
 
 
 
◆ showStack()
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder.showStack | ( |  | binLattices, |  
          |  |  |  | plane = 'xz', |  
          |  |  |  | scl = None, |  
          |  |  |  | maxWidth = 600, |  
          |  |  |  | maxHeight = 600 ) |  | static | 
 
Display all slices of the site lattice interactively
Args:
    binLattices:
        Binary lattice or collection of binary lattices
        Viewing plane, e.g. "xy"
    plane (str):
        Plane to display: (xy, yz, xz)
    scl (int):
        Scale pixels by this amount
    maxWidth (int):
        Maximum width of image
    maxHeight (int):
        Maximum height of image
 
 
 
◆ sphereStructElem()
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder.sphereStructElem | ( |  | r | ) |  |  | static | 
 
Spherical structuring element
Args:
    r (int):
        radius
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Structuring element
 
 
 
◆ spoke()
      
        
          | jLM.RegionBuilder.RegionBuilder.spoke | ( |  | self, | 
        
          |  |  |  | x0, | 
        
          |  |  |  | length, | 
        
          |  |  |  | spoke_radius, | 
        
          |  |  |  | r, | 
        
          |  |  |  | phi, | 
        
          |  |  |  | theta ) | 
      
 
Construct spoke
For a sphere of radius, :math:`r`, centered on :math:`x_0`, a spoke is
a cylinder centered on the surface of the sphere at the polar
coordinates :math:`(r,\phi,\theta)`, rotated to be normal to the sphere
surface.
Args:
    x0 (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        Center of sphere
    length (float):
        Length of spoke. Spoke will protrude 0.5*length from the inside and outside of the sphere
    spoke_radius (float):
        Radius (thickness) of the spoke
    r (float):
        Radius of sphere
    phi (float):
        Azimuthal position of spoke with respect to the sphere center
    theta (float):
        Polar position of spoke with respect to the sphere center
Returns:
    :py:class:`~numpy.ndarray(shape=(nx,ny,nz), dtype=bool)`: 
        Lattice mask
 
 
 
◆ transformGrid()
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder.transformGrid | ( |  | xs, |  
          |  |  |  | x0, |  
          |  |  |  | alpha, |  
          |  |  |  | beta, |  
          |  |  |  | gamma ) |  | static | 
 
Compute the translation/rotation of an index grid
Args:
    xs (:py:class:`~numpy.ndarray(shape=(3,nx,ny,nz), dtype=int)`): 
        Index grid of dimensions [3, nx, ny, nz]
    x0 (:py:class:`~numpy.ndarray(shape=(3,), dtype=float)`):
        Center of rotation
    alpha (float):
        Euler Z-rotation (radians)
    beta (float):
        Euler X-rotation (radians)
    gamma (float):
        Euler Z'-rotation (radians)
Returns:
    :py:class:`~numpy.ndarray(shape=(3,nx,ny,nz), dtype=int)`: 
        Transformed grid
 
 
 
◆ _net
  
  | 
        
          | jLM.RegionBuilder.RegionBuilder._net = net |  | protected | 
 
 
◆ center
      
        
          | jLM.RegionBuilder.RegionBuilder.center = np.array( [self.nx//2, self.ny//2, self.nz//2]) | 
      
 
 
◆ nx
      
        
          | jLM.RegionBuilder.RegionBuilder.nx | 
      
 
 
◆ ny
      
        
          | jLM.RegionBuilder.RegionBuilder.ny | 
      
 
 
◆ nz
      
        
          | jLM.RegionBuilder.RegionBuilder.nz = net.siteLattice.shape | 
      
 
 
◆ origin
      
        
          | jLM.RegionBuilder.RegionBuilder.origin = np.zeros( 3 ) | 
      
 
 
◆ shape
      
        
          | jLM.RegionBuilder.RegionBuilder.shape = net.siteLattice.shape | 
      
 
 
◆ xs
      
        
          | jLM.RegionBuilder.RegionBuilder.xs = np.mgrid[0:self.nx, 0:self.ny, 0:self.nz] | 
      
 
 
The documentation for this class was generated from the following file: