Skip to content

conversion

Turn mesh neuron into skeleton.

This function is a thin-wrapper for skeletor. It uses sensible defaults for neurons but if you want to fine-tune your skeletons you should look into using skeletor directly.

PARAMETER DESCRIPTION
x
    Mesh(es) to skeletonize. Note that the quality of the results
    very much depends on the mesh, so it might be worth doing some
    pre-processing (see below).

TYPE: MeshNeuron | trimesh.Trimesh

method
    Method to use for skeletonization:
     - "wavefront": fast but noisier, skeletons will be ~centered
       within the neuron
     - "teasar": slower but smoother, skeletons follow the
       surface of the mesh, requires the `inv_dist` parameter to be
       set
    "wavefront" also produces radii, "teasar" doesn't.

TYPE: 'wavefront' | 'teasar' DEFAULT: 'wavefront'

fix_mesh
    Whether to try to fix some common issues in the mesh before
    skeletonization. Note that this might compromise the
    vertex-to-node-ID mapping.

TYPE: bool DEFAULT: False

shave
    Whether to "shave" the resulting skeleton to reduce bristles
    on the backbone.

TYPE: bool DEFAULT: True

heal
    Whether to heal the resulting skeleton if it is fragmented.
    For more control over the stitching set `heal=False` and use
    [`navis.heal_skeleton`][] directly. Note that this
    can be fairly costly if the mesh as many tiny fragments.

TYPE: bool | "LEAFS" | "ALL" DEFAULT: False

connectors
    Whether to carry over existing connector tables. This will
    attach connectors by first snapping them to the closest mesh
    vertex and then to the skeleton node corresponding to that
    vertex.

TYPE: bool DEFAULT: False

inv_dist
    Only required for method "teasar": invalidation distance for
    the traversal. Smaller `inv_dist` captures smaller features
    but is slower and more noisy, and vice versa. A good starting
    value is around 2-5 microns. Can be a unit string - e.g.
    "5 microns" - if your neuron has its units set.

TYPE: int | float | str DEFAULT: None

**kwargs
    Additional keyword arguments are passed through to the respective
    function in `skeletor` - i.e. `by_wavefront` or `by_teasar`.

DEFAULT: {}

RETURNS DESCRIPTION
skeleton

Has a .vertex_map attribute that maps each vertex in the input mesh to a skeleton node ID.

TYPE: navis.TreeNeuron

See Also

navis.drop_fluff Use this if your mesh has lots of tiny free floating bits to reduce noise and speed up skeletonization.

Examples:

>>> import navis
>>> # Get a mesh neuron
>>> n = navis.example_neurons(1, kind='mesh')
>>> # Convert to skeleton
>>> sk = navis.conversion.mesh2skeleton(n)
>>> # Mesh vertex indices to node IDs map
>>> sk.vertex_map
array([938, 990, 990, ...,  39, 234, 234])
Source code in navis/conversion/converters.py
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
@utils.map_neuronlist(desc='Skeletonizing', allow_parallel=True)
def mesh2skeleton(x: 'core.MeshNeuron',
                  method: str = 'wavefront',
                  fix_mesh: bool = False,
                  shave: bool = True,
                  heal: bool = False,
                  connectors: bool = False,
                  inv_dist: Union[int, float] = None,
                  **kwargs):
    """Turn mesh neuron into skeleton.

    This function is a thin-wrapper for `skeletor`. It uses sensible defaults
    for neurons but if you want to fine-tune your skeletons you should look
    into using `skeletor` directly.

    Parameters
    ----------
    x :         MeshNeuron | trimesh.Trimesh
                Mesh(es) to skeletonize. Note that the quality of the results
                very much depends on the mesh, so it might be worth doing some
                pre-processing (see below).
    method :    'wavefront' | 'teasar'
                Method to use for skeletonization:
                 - "wavefront": fast but noisier, skeletons will be ~centered
                   within the neuron
                 - "teasar": slower but smoother, skeletons follow the
                   surface of the mesh, requires the `inv_dist` parameter to be
                   set
                "wavefront" also produces radii, "teasar" doesn't.
    fix_mesh :  bool
                Whether to try to fix some common issues in the mesh before
                skeletonization. Note that this might compromise the
                vertex-to-node-ID mapping.
    shave :     bool
                Whether to "shave" the resulting skeleton to reduce bristles
                on the backbone.
    heal :      bool | "LEAFS" | "ALL"
                Whether to heal the resulting skeleton if it is fragmented.
                For more control over the stitching set `heal=False` and use
                [`navis.heal_skeleton`][] directly. Note that this
                can be fairly costly if the mesh as many tiny fragments.
    connectors : bool
                Whether to carry over existing connector tables. This will
                attach connectors by first snapping them to the closest mesh
                vertex and then to the skeleton node corresponding to that
                vertex.
    inv_dist :  int | float | str
                Only required for method "teasar": invalidation distance for
                the traversal. Smaller `inv_dist` captures smaller features
                but is slower and more noisy, and vice versa. A good starting
                value is around 2-5 microns. Can be a unit string - e.g.
                "5 microns" - if your neuron has its units set.
    **kwargs
                Additional keyword arguments are passed through to the respective
                function in `skeletor` - i.e. `by_wavefront` or `by_teasar`.

    Returns
    -------
    skeleton :  navis.TreeNeuron
                Has a `.vertex_map` attribute that maps each vertex in the
                input mesh to a skeleton node ID.

    See Also
    --------
    [`navis.drop_fluff`][]
                Use this if your mesh has lots of tiny free floating bits to
                reduce noise and speed up skeletonization.

    Examples
    --------
    >>> import navis
    >>> # Get a mesh neuron
    >>> n = navis.example_neurons(1, kind='mesh')
    >>> # Convert to skeleton
    >>> sk = navis.conversion.mesh2skeleton(n)
    >>> # Mesh vertex indices to node IDs map
    >>> sk.vertex_map                                           # doctest: +SKIP
    array([938, 990, 990, ...,  39, 234, 234])

    """
    utils.eval_param(x, name='x', allowed_types=(core.MeshNeuron, tm.Trimesh))
    utils.eval_param(method, name='method', allowed_values=('wavefront', 'teasar'))

    if method == 'teasar' and inv_dist is None:
        raise ValueError('Must set `inv_dist` parameter when using method '
                         '"teasar". A good starting value is around 2-5 microns.')

    props = {'soma': None}
    if isinstance(x, core.MeshNeuron):
        props.update({'id': x.id, 'name': x.name, 'units': x.units})
        if x.has_soma_pos:
            props['soma_pos'] = x.soma_pos

        if not isinstance(inv_dist, type(None)):
            inv_dist = x.map_units(inv_dist)

        mesh = x.trimesh
    else:
        mesh = x
        x = core.MeshNeuron(x)

    if fix_mesh:
        mesh = sk.pre.fix_mesh(mesh, remove_disconnected=False)

    kwargs['progress'] = False
    if method == 'wavefront':
        skeleton = sk.skeletonize.by_wavefront(mesh, **kwargs)
    elif method == 'teasar':
        skeleton = sk.skeletonize.by_teasar(x, inv_dist=inv_dist, **kwargs)

    props['vertex_map'] = skeleton.mesh_map

    s = core.TreeNeuron(skeleton.swc, **props)

    if s.has_soma:
        s.reroot(s.soma, inplace=True)

    if heal:
        _ = morpho.heal_skeleton(s, inplace=True, method='ALL')

    if shave:
        # Find single node bristles
        leafs = s.leafs.node_id.values

        # Make sure we keep the soma
        if s.has_soma:
            leafs = leafs[~np.isin(leafs, s.soma)]

        bp = s.branch_points.node_id.values
        bristles = s.nodes[s.nodes.node_id.isin(leafs)
                           & s.nodes.parent_id.isin(bp)]

        # Subset neuron
        keep = s.nodes[~s.nodes.node_id.isin(bristles.node_id)].node_id.values
        s = morpho.subset_neuron(s, keep, inplace=True)

        # Fix vertex map
        for b, p in zip(bristles.node_id.values, bristles.parent_id.values):
            s.vertex_map[s.vertex_map == b] = p

    # In particular with method wavefront, some nodes (mostly leafs) can have
    # a radius of 0. We will fix this here by giving them 1/2 the radius of
    # their parent nodes'
    to_fix = (s.nodes.radius == 0) & (s.nodes.parent_id >= 0)
    if any(to_fix):
        radii = s.nodes.set_index('node_id').radius
        new_radii = radii.loc[s.nodes.loc[to_fix].parent_id].values / 2
        s.nodes.loc[to_fix, 'radius'] = new_radii

    # Last but not least: map connectors
    if connectors and x.has_connectors:
        cn_table = x.connectors.copy()

        # A connector/id column is currently required for skeletons but not
        # meshes
        if not any(np.isin(('id', 'connector_id'), cn_table.columns)):
            cn_table.insert(0, 'connector_id', np.arange(len(cn_table)))

        cn_table['node_id'] = x.snap(cn_table[['x', 'y', 'z']].values)[0]
        node_map = dict(zip(np.arange(len(s.vertex_map)), s.vertex_map))
        cn_table['node_id'] = cn_table.node_id.map(node_map)
        s.connectors = cn_table

    return s

Turn neuron into voxels.

PARAMETER DESCRIPTION
x
        Neuron(s) to voxelize. Uses the neurons' nodes, vertices and
        points, respectively.

TYPE: TreeNeuron | MeshNeuron | Dotprops

pitch
        Side length(s) voxels. Can be isometric (float) or an
        iterable of dimensions in (x, y, z).

TYPE: float | iterable thereof

bounds
        Boundaries [in units of `x`] for the voxel grid. If not
        provided, will use `x.bbox`.

TYPE: (3, 2) or (2, 3) array DEFAULT: None

counts
        If True, voxel grid will have point counts for values
        instead of just True/False.

TYPE: bool DEFAULT: False

vectors
        If True, will also attach a vector field as `.vectors`
        property.

TYPE: bool DEFAULT: False

alphas
        If True, will also return a grid with alpha values as
        `.alpha` property.

TYPE: bool DEFAULT: False

smooth
        If non-zero, will apply a Gaussian filter with `smooth`
        as `sigma`.

TYPE: int DEFAULT: 0

RETURNS DESCRIPTION
VoxelNeuron

Has the voxel grid as .grid and (optionally) .vectors and .alphas properties. .grid data type depends on settings: - default = bool (i.e. True/False) - if counts=True = integer - if smooth=True = float Empty voxels will have vector (0, 0, 0) and alpha 0. Also note that data tables (e.g. connectors) are not carried over from the input neuron.

Examples:

>>> import navis
>>> # Get a skeleton
>>> n = navis.example_neurons(1)
>>> # Convert to voxel neuron
>>> vx = navis.conversion.neuron2voxels(n, pitch='5 microns')
Source code in navis/conversion/converters.py
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
@utils.map_neuronlist(desc='Voxelizing', allow_parallel=True)
def neuron2voxels(x: 'core.BaseNeuron',
                  pitch: Union[list, tuple, float],
                  bounds: Optional[list] = None,
                  counts: bool = False,
                  vectors: bool = False,
                  alphas: bool = False,
                  smooth: int = 0) -> 'core.VoxelNeuron':
    """Turn neuron into voxels.

    Parameters
    ----------
    x :             TreeNeuron | MeshNeuron | Dotprops
                    Neuron(s) to voxelize. Uses the neurons' nodes, vertices and
                    points, respectively.
    pitch :         float | iterable thereof
                    Side length(s) voxels. Can be isometric (float) or an
                    iterable of dimensions in (x, y, z).
    bounds :        (3, 2)  or (2, 3) array, optional
                    Boundaries [in units of `x`] for the voxel grid. If not
                    provided, will use `x.bbox`.
    counts :        bool
                    If True, voxel grid will have point counts for values
                    instead of just True/False.
    vectors :       bool
                    If True, will also attach a vector field as `.vectors`
                    property.
    alphas :        bool
                    If True, will also return a grid with alpha values as
                    `.alpha` property.
    smooth :        int
                    If non-zero, will apply a Gaussian filter with `smooth`
                    as `sigma`.

    Returns
    -------
    VoxelNeuron
                    Has the voxel grid as `.grid` and (optionally) `.vectors`
                    and `.alphas` properties. `.grid` data type depends
                    on settings:
                     - default = bool (i.e. True/False)
                     - if `counts=True` = integer
                     - if `smooth=True` = float
                    Empty voxels will have vector (0, 0, 0) and alpha 0. Also
                    note that data tables (e.g. `connectors`) are not carried
                    over from the input neuron.

    Examples
    --------
    >>> import navis
    >>> # Get a skeleton
    >>> n = navis.example_neurons(1)
    >>> # Convert to voxel neuron
    >>> vx = navis.conversion.neuron2voxels(n, pitch='5 microns')

    """
    if not utils.is_iterable(pitch):
        # Map units (non-str are just passed through)
        pitch = x.map_units(pitch, on_error='raise')
        if not isinstance(pitch, Number):
            raise TypeError('Expected `pitch` to be a number (or list thereof)'
                            f', got {type(pitch)}')
        pitch = [pitch] * 3
    elif len(pitch) != 3:
        raise ValueError('`pitch` must be single number or a list of three')
    else:
        pitch = np.array([x.map_units(p, on_error='raise') for p in pitch])

    # Convert to voxel indices
    ix, _ = _make_voxels(x=x, pitch=pitch, strip=False)

    if isinstance(bounds, type(None)):
        bounds = x.bbox
    else:
        bounds = np.asarray(bounds)

    if bounds.shape == (2, 3):
        bounds = bounds.T

    # Shape of grid
    dim = np.ceil(bounds[:, 1] / pitch) - np.floor(bounds[:, 0] / pitch)
    shape = np.ceil(dim).astype(int) + 1

    # Get unique voxels
    if not counts:
        vxl = np.unique(ix, axis=0)
    else:
        vxl, cnt = np.unique(ix, axis=0, return_counts=True)

    # Substract lower bounds
    offset = (bounds[:, 0] / pitch)
    vxl = vxl - offset.round().astype(int)
    ix = ix - offset.round().astype(int)

    # Drop voxels outside the defined bounds
    vxl = vxl[vxl.min(axis=1) >= 0]
    vxl = vxl[np.all(vxl < shape, axis=1)]

    # Generate grid
    grid = np.zeros(shape=shape, dtype=bool)

    # Populate grid
    if not counts:
        grid[vxl[:, 0], vxl[:, 1], vxl[:, 2]] = True
    else:
        grid = grid.astype(int)
        grid[vxl[:, 0], vxl[:, 1], vxl[:, 2]] = cnt

    # Apply Gaussian filter
    if smooth:
        grid = gaussian_filter(grid.astype(np.float32), sigma=smooth)

    # Generate neuron
    units = [f'{p * u} {x.units.units}' for p, u in zip(utils.make_iterable(pitch),
                                                        x.units_xyz.magnitude)]
    offset = offset * pitch * x.units_xyz.magnitude
    n = core.VoxelNeuron(grid, id=x.id, name=x.name, units=units, offset=offset)

    # If no vectors required, we can just return now
    if not vectors and not alphas:
        return n

    if isinstance(x, core.TreeNeuron):
        pts = x.nodes[['x', 'y', 'z']].values
    elif isinstance(x, core.Dotprops):
        pts = x.points
    elif isinstance(x, core.MeshNeuron):
        pts = np.array(x.vertices)

    # Generate an empty vector field
    vects = np.zeros((grid.shape[0], grid.shape[1], grid.shape[2], 3),
                     dtype=np.float32)
    alph = np.zeros(grid.shape, dtype=np.float32)

    # Get unique voxels
    uni, inv = np.unique(ix, axis=0, return_inverse=True)

    # Go over each voxel
    for i in range(len(uni)):
        # Get points in this voxel
        pt = pts[inv == i]

        # Reshape
        pt = pt.reshape(1, -1, 3)

        # Generate centers for each cloud of k nearest neighbors
        centers = np.mean(pt, axis=1)

        # Generate vector from center
        cpt = pt - centers.reshape((pt.shape[0], 1, 3))

        # Get inertia (N, 3, 3)
        inertia = cpt.transpose((0, 2, 1)) @ cpt

        # Extract vector and alpha
        u, s, vh = np.linalg.svd(inertia)
        vect = vh[:, 0, :]

        # No alpha if only one point
        if pt.shape[1] > 1:
            alpha = (s[:, 0] - s[:, 1]) / np.sum(s, axis=1)
        else:
            alpha = [0]

        vects[uni[i][0], uni[i][1], uni[i][2]] = vect.flatten()
        alph[uni[i][0], uni[i][1], uni[i][2]] = alpha[0]

    if vectors:
        n.vectors = vects
    if alpha:
        n.alphas = alpha

    return n

Turn points into skeleton.

This function works by: 1. Compute the k nearest neighbors for each point 2. Generate a graph from the nearest-neighbor edges 3. Extract a minimum-spanning tree (MST) from the graph 4. Process the MST into a skeleton

PARAMETER DESCRIPTION
x
    Points to skeletonize.

TYPE: (N, 3) array | Dotprops

k
    Number of nearest neighbors to consider. Too low values of `k`
    can lead to disconnected skeletons.

TYPE: int DEFAULT: 10

max_dist
    Edges longer than this will be ignored. This can lead to a
    fragmented (i.e. multi-root) skeleton!

TYPE: float DEFAULT: None

RETURNS DESCRIPTION
skeleton

TYPE: navis.TreeNeuron

Examples:

>>> import navis
>>> # Get a mesh neuron
>>> n = navis.example_neurons(1)
>>> # Get the points
>>> pts = n.nodes[['x', 'y', 'z']].values
>>> # Convert points back into skeleton
>>> sk = navis.conversion.points2skeleton(pts)
Source code in navis/conversion/converters.py
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
@utils.map_neuronlist(desc='Skeletonizing', allow_parallel=True)
def points2skeleton(x: Union['core.Dotprops', np.ndarray],
                    k: int = 10,
                    max_dist: Optional[float] = None):
    """Turn points into skeleton.

    This function works by:
     1. Compute the `k` nearest neighbors for each point
     2. Generate a graph from the nearest-neighbor edges
     3. Extract a minimum-spanning tree (MST) from the graph
     4. Process the MST into a skeleton

    Parameters
    ----------
    x :         (N, 3) array | Dotprops
                Points to skeletonize.
    k :         int
                Number of nearest neighbors to consider. Too low values of `k`
                can lead to disconnected skeletons.
    max_dist :  float, optional
                Edges longer than this will be ignored. This can lead to a
                fragmented (i.e. multi-root) skeleton!

    Returns
    -------
    skeleton :  navis.TreeNeuron

    Examples
    --------
    >>> import navis
    >>> # Get a mesh neuron
    >>> n = navis.example_neurons(1)
    >>> # Get the points
    >>> pts = n.nodes[['x', 'y', 'z']].values
    >>> # Convert points back into skeleton
    >>> sk = navis.conversion.points2skeleton(pts)

    """
    utils.eval_param(x, name='x', allowed_types=(core.Dotprops, np.ndarray))

    if isinstance(x, core.Dotprops):
        pts = x.points
    else:
        if (x.ndim != 2) and (x.shape[1] != 3):
            raise ValueError(f'Points must be shape (N, 3), got {x.shape}')
        pts = x

    # Get the list of nearest neighbours
    tree = core.dotprop.KDTree(pts)

    defaults = {}
    if max_dist is not None:
        # We have to avoid passing `None` because scipy's KDTree does not like
        # that (pykdtree does not care)
        defaults['distance_upper_bound'] = max_dist
    dists, NN = tree.query(pts, k=k + 1, **defaults)

    # Drop self-hits
    dists, NN = dists[:, 1:], NN[:, 1:]

    # Turn into edges
    edges = []
    ix1 = np.arange(len(dists))
    for i in range(k):
        ix2 = NN[:, i]
        le = dists[:, i]
        # If a max dist was set we have to remove NN that have dist np.inf
        if max_dist is None:
            edges += list(zip(ix1, ix2, le))
        else:
            not_inf = le != np.inf
            edges += list(zip(ix1[not_inf], ix2[not_inf], le[not_inf]))

    # Generate graph
    G = nx.Graph()
    G.add_nodes_from(ix1)
    G.add_weighted_edges_from(edges)

    # Extract minimum spanning tree
    G_mst = nx.minimum_spanning_tree(G)

    # Add the coordinates as node properties
    nx.set_node_attributes(G_mst, dict(zip(G.nodes, pts[:, 0])), name='x')
    nx.set_node_attributes(G_mst, dict(zip(G.nodes, pts[:, 1])), name='y')
    nx.set_node_attributes(G_mst, dict(zip(G.nodes, pts[:, 2])), name='z')

    return graph.nx2neuron(G_mst)

Convert TreeNeuron to MeshNeuron.

Uses the radius to convert skeleton to 3D tube mesh. Missing radii are treated as zeros.

PARAMETER DESCRIPTION
x
        Neuron to convert.

TYPE: TreeNeuron | NeuronList

tube_points
        Number of points making up the circle of the cross-section
        of the tube.

TYPE: int DEFAULT: 8

radius_scale_factor
        Factor to scale radii by.

TYPE: float DEFAULT: 1

use_normals
        If True will rotate tube along its curvature.

TYPE: bool DEFAULT: True

warn_missing_radii
        Whether to warn if radii are missing or <= 0.

TYPE: bool DEFAULT: True

RETURNS DESCRIPTION
TreeNeuron

Data tables (e.g. connectors) are not carried over from the input neuron.

Examples:

>>> import navis
>>> # Get a skeleton
>>> n = navis.example_neurons(1)
>>> # Convert to mesh neuron
>>> m = navis.conversion.tree2meshneuron(n)
Source code in navis/conversion/converters.py
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
@utils.map_neuronlist(desc='Converting', allow_parallel=True)
def tree2meshneuron(x: 'core.TreeNeuron',
                    tube_points: int = 8,
                    radius_scale_factor: float = 1,
                    use_normals: bool = True,
                    warn_missing_radii: bool = True
                    ) -> 'core.MeshNeuron':
    """Convert TreeNeuron to MeshNeuron.

    Uses the `radius` to convert skeleton to 3D tube mesh. Missing radii are
    treated as zeros.

    Parameters
    ----------
    x :             TreeNeuron | NeuronList
                    Neuron to convert.
    tube_points :   int
                    Number of points making up the circle of the cross-section
                    of the tube.
    radius_scale_factor : float
                    Factor to scale radii by.
    use_normals :   bool
                    If True will rotate tube along its curvature.
    warn_missing_radii : bool
                    Whether to warn if radii are missing or <= 0.

    Returns
    -------
    TreeNeuron
                    Data tables (e.g. `connectors`) are not carried over from
                    the input neuron.

    Examples
    --------
    >>> import navis
    >>> # Get a skeleton
    >>> n = navis.example_neurons(1)
    >>> # Convert to mesh neuron
    >>> m = navis.conversion.tree2meshneuron(n)

    """
    # Delay to avoid circular imports
    from ..plotting.plot_utils import make_tube

    if not isinstance(x, core.TreeNeuron):
        raise TypeError(f'Expected TreeNeuron, got "{type(x)}"')

    # Map segments of node IDs to segments of node indices
    id2ix = dict(zip(x.nodes.node_id, np.arange(len(x.nodes))))
    segments = [np.array([id2ix[n] for n in seg]) for seg in x.segments]

    # Note that we are treating missing radii as "0"
    radii_map = x.nodes.radius.fillna(0).values
    if warn_missing_radii and (radii_map <= 0).any():
        logger.warning('At least some radii are missing or <= 0. Mesh may look funny.')

    # Map radii onto segments
    radii = [radii_map[seg] * radius_scale_factor for seg in segments]
    co_map = x.nodes[['x', 'y', 'z']].values
    seg_points = [co_map[seg] for seg in segments]

    vertices, faces = make_tube(seg_points,
                                radii=radii,
                                tube_points=tube_points,
                                use_normals=use_normals)

    # Note: the `process=False` is necessary to not break correspondence
    # by e.g. merging duplicate vertices
    m = core.MeshNeuron({'vertices': vertices, 'faces': faces},
                        units=x.units, name=x.name, id=x.id, process=False)


    # For each vertex, track the original node: the first `tube_points` vertices
    # correspond to the first node of the first segment and so on.
    m.vertex_map = np.concatenate([np.repeat(seg, tube_points) for seg in segments])

    return m

Generate mesh from voxels using marching cubes.

PARAMETER DESCRIPTION
voxels
        Object to voxelize. Can be a VoxelNeuron or an (N, 3) array
        of x, y, z voxel coordinates.

TYPE: VoxelNeuron | (N, 3) np.array

spacing
        (3, ) array with x, y, z voxel size. If `auto` and input is
        a `VoxelNeuron` we will use the neuron's `.units` property,
        else spacing will be `(1, 1, 1)`.

TYPE: np.array DEFAULT: 'auto'

step_size
        Step size for marching cube algorithm.
        Higher values = faster but coarser.

TYPE: int DEFAULT: 1

chunk_size
        Whether to process voxels in chunks to keep memory footprint
        low:
          - "auto" will set chunk size automatically based on size
            of input
          - use `int` to set chunk size - smaller chunk mean lower
            memory consumption but longer run time - 200 (i.e.
            chunks of 200x200x200 voxels) appears to be a good value
          - set to `0` to force processing in one go

TYPE: "auto" | int DEFAULT: 'auto'

For

pad_chunks
        If True, will pad each chunk. This helps making meshes
        watertight but may introduce internal faces when merging
        mesh fragments.

TYPE: bool DEFAULT: True

merge_fragments
        If True, will attempt to merge fragments at the chunk
        boundaries.

TYPE: bool DEFAULT: True

progress
        Whether to show a progress bar.

TYPE: bool DEFAULT: True

RETURNS DESCRIPTION
mesh

Returns a trimesh or MeshNeuron depending on the input. Data tables (e.g. connectors) are not carried over from input neuron.

TYPE: trimesh.Trimesh | MeshNeuron

Source code in navis/conversion/meshing.py
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
def voxels2mesh(vox: Union['core.VoxelNeuron', np.ndarray],
                spacing: Union[Literal['auto'], np.ndarray] = 'auto',
                step_size: int = 1,
                chunk_size: Optional[int] = 'auto',
                pad_chunks: bool = True,
                merge_fragments: bool = True,
                progress: bool = True) -> Union[tm.Trimesh, 'core.MeshNeuron']:
    """Generate mesh from voxels using marching cubes.

    Parameters
    ----------
    voxels :        VoxelNeuron | (N, 3) np.array
                    Object to voxelize. Can be a VoxelNeuron or an (N, 3) array
                    of x, y, z voxel coordinates.
    spacing :       np.array
                    (3, ) array with x, y, z voxel size. If `auto` and input is
                    a `VoxelNeuron` we will use the neuron's `.units` property,
                    else spacing will be `(1, 1, 1)`.
    step_size :     int, optional
                    Step size for marching cube algorithm.
                    Higher values = faster but coarser.
    chunk_size :    "auto" | int, optional
                    Whether to process voxels in chunks to keep memory footprint
                    low:
                      - "auto" will set chunk size automatically based on size
                        of input
                      - use `int` to set chunk size - smaller chunk mean lower
                        memory consumption but longer run time - 200 (i.e.
                        chunks of 200x200x200 voxels) appears to be a good value
                      - set to `0` to force processing in one go

    For `chunk_size != 0`:

    pad_chunks :    bool
                    If True, will pad each chunk. This helps making meshes
                    watertight but may introduce internal faces when merging
                    mesh fragments.
    merge_fragments :  bool
                    If True, will attempt to merge fragments at the chunk
                    boundaries.
    progress :      bool
                    Whether to show a progress bar.

    Returns
    -------
    mesh :          trimesh.Trimesh | MeshNeuron
                    Returns a trimesh or MeshNeuron depending on the input.
                    Data tables (e.g. `connectors`) are not carried over from
                    input neuron.

    """
    if not skimage:
        raise ModuleNotFoundError(
            'Meshing requires `skimage`:\n '
            'pip3 install scikit-image'
            )

    utils.eval_param(vox, 'vox', allowed_types=(core.VoxelNeuron, np.ndarray))

    if spacing == 'auto':
        if not isinstance(vox, core.VoxelNeuron):
            spacing = [1, 1, 1]
        else:
            spacing = vox.units_xyz.magnitude

    if isinstance(vox, core.VoxelNeuron):
        voxels = vox.voxels
    else:
        voxels = vox

    if voxels.ndim != 2 or voxels.shape[1] != 3:
        raise ValueError(f'Voxels must be shape (N, 3), got {voxels.shape}')

    if chunk_size == 'auto':
        if len(voxels) > 1e6:
            chunk_size = 200
        else:
            chunk_size = 0

    if not chunk_size:
        mesh = _mesh_from_voxels_single(voxels=voxels,
                                        spacing=spacing,
                                        step_size=step_size)
    else:
        mesh = _mesh_from_voxels_chunked(voxels=voxels,
                                         spacing=spacing,
                                         step_size=step_size,
                                         chunk_size=chunk_size,
                                         pad_chunks=pad_chunks,
                                         merge_fragments=merge_fragments,
                                         progress=progress)

    if isinstance(vox, core.VoxelNeuron):
        mesh.vertices += vox.offset
        mesh = core.MeshNeuron(mesh, units=f'1 {vox.units.units}', id=vox.id)

    return mesh