fiber module

Defines the Fiber class and helper functions for building fiber models.

This module provides functionality for building and simulating both 1D and 3D fiber models in the NEURON environment.

pyfibers.fiber.build_fiber(fiber_model, diameter, length=None, n_sections=None, n_nodes=None, enforce_odd_nodecount=True, **kwargs)

Generate a 1D (straight) fiber model in NEURON.

This function creates a model fiber as an instance of the Fiber class using the specific subclass specified from the FiberModel enumerator. with user-specified diameter and length (from one of: number of sections, number of nodes, or length in microns). Additional keyword arguments are forwarded to the fiber model class constructor.

By default, the first section of the fiber is located at the origin (0, 0, 0), and the fiber extends along the z-axis in the positive direction. To change the fiber’s location, the method Fiber.set_xyz() can be used to translate the fiber along the x, y, or z axes. To create fibers along a custom path in 3D space, use build_fiber_3d() instead.

Parameters:
  • fiber_model (FiberModel) – A FiberModel enumerator specifying the type of fiber to instantiate.

  • diameter (float) – The fiber diameter in micrometers (µm).

  • length (float) – The total length of the fiber in micrometers (µm), if defining by length.

  • n_sections (int) – The total number of sections for discretizing the fiber, if defining by sections.

  • n_nodes (int) – The total number of nodes along the fiber, if defining by nodes.

  • enforce_odd_nodecount (bool) – If True, ensure that the number of nodes is odd.

  • kwargs – Additional arguments forwarded to the underlying fiber model class.

Raises:

ValueError – If more than one among length, n_sections, or n_nodes is specified.

Return type:

Fiber

Returns:

A Fiber class instance.

Example Usage

from PyFibers import build_fiber, FiberModel

fiber = build_fiber(fiber_model=FiberModel.MRG_INTERPOLATION, diameter=10, n_nodes=25)
pyfibers.fiber.build_fiber_3d(fiber_model, diameter, path_coordinates, shift=0, shift_ratio=None, center_shift=False, **kwargs)

Generate a 3D fiber model in NEURON based on a specified path.

This function calculates the fiber’s length from the user-supplied path_coordinates and uses it internally to instantiate a 3D fiber model. The coordinates are a 2D numpy array of shape (number_of_points, 3), where each row represents a point in 3D space (x, y, z).

The fiber model will be created by repeating sections along the path until no more nodes can be added without exceeding the path length. By default, the center of the first section is placed at the origin (0, 0, 0), and the fiber extends along the path.

Parameters:
  • fiber_model (FiberModel) – A FiberModel enumerator specifying the type of fiber to instantiate.

  • diameter (float) – The fiber diameter in micrometers (µm).

  • path_coordinates (ndarray) – A numpy array of shape (N, 3) specifying the 3D coordinates (x, y, z) of the fiber path.

  • shift (float) – A shift in microns to apply to the fiber coordinates.

  • shift_ratio (float) – Ratio of the internodal length to shift the fiber coordinates.

  • center_shift (bool) – If True, center the fiber before applying the shift.

  • kwargs – Additional arguments forwarded to the underlying fiber model class.

Raises:

ValueError – If path_coordinates is not provided, or if n_sections, n_nodes, or length is specified (these are invalid in 3D mode).

Return type:

Fiber

Returns:

A fully instantiated 3D fiber model Fiber instance.

Example:

import numpy as np
from PyFibers import build_fiber_3d, FiberModel

coords = np.array(
    [
        [0.0, 0.0, 0.0],
        [0.0, 1.0, 3.0],
        [0.0, 2.0, 7.0],
        # ...
    ]
)
fiber = build_fiber_3d(
    fiber_model=FiberModel.MRG_INTERPOLATION, diameter=10, path_coordinates=coords
)
print(fiber)
class pyfibers.fiber.Fiber(fiber_model, diameter, temperature=37, passive_end_nodes=True, is_3d=False)

Base class for model fibers.

The Fiber class provides functionality for constructing, configuring, and simulating fiber models. It encapsulates key methods for:

  • Generating sections specified by a fiber model subclass

  • Recording membrane voltage, current, and gating variables

  • Calculating extracellular potentials and single fiber action potentials

  • Measuring conduction velocity along the fiber

  • Handling 3D or 1D fiber geometry

Initialize the Fiber class.

Parameters:
  • fiber_model (FiberModel) – The enumerator attribute (from FiberModel) representing the type of fiber model.

  • diameter (float) – The diameter of the fiber (µm).

  • temperature (float) – The temperature at which the fiber will be simulated, in Celsius.

  • passive_end_nodes (int | bool) – If True, automatically assign passive properties to the end nodes. Can also be an integer specifying how many passive end nodes to include at each end.

  • is_3d (bool) – If True, fiber coordinates are treated as 3D. Usually set automatically by build_fiber_3d().

Variables:
  • v_rest – The resting membrane potential of the fiber (in mV).

  • gating_variables – A dictionary mapping gating variable names to their corresponding NEURON attribute names.

  • myelinatedTrue if myelinated, False if unmyelinated.

Variables:
  • diameter – The diameter of the fiber in micrometers.

  • fiber_model – The FiberModel attribute name.

  • temperature – The temperature at which the fiber will be simulated [C].

  • passive_end_nodes – The number of passive end nodes included at each end.

  • delta_z – The center-to-center internodal length of the fiber.

Variables:
  • nodecount – The number of nodes in the fiber.

  • length – The total length of the fiber in micrometers (end-to-end).

  • sections – A list of NEURON Section objects representing the fiber.

  • nodes – A list of NEURON Section objects representing the fiber nodes only (subset of Fiber.sections).

  • coordinates – A numpy array of 3D coordinates for the center of each section along the fiber.

  • longitudinal_coordinates – A numpy array of 1D (arc-length) coordinates of the center of each section along the fiber.

  • path – A nd_line object representing the 3D path of the fiber.

Variables:
  • apc – A list of NEURON APCount objects for detecting action potentials (from :meth:Fiber.apcounts()).

  • vm – A list of NEURON Vector objects recording membrane voltage (from :meth:Fiber.record_vm()).

  • im – A list of NEURON Vector objects recording membrane current (from :meth:Fiber.record_im()).

  • vext – A list of NEURON Vector objects recording extracellular potential (from :meth:Fiber.record_vext()).

  • time – A NEURON Vector recording the simulation time.

  • gating – A dictionary mapping gating variable names to lists of recorded NEURON Vector objects (from :meth:Fiber.record_gating()).

Variables:
  • nc – A NEURON NetCon object for intrinsic activity.

  • syn – A NEURON ExpSyn object for intrinsic activity.

  • stim – A NEURON NetStim object for intrinsic activity.

Variables:

potentials – A numpy array of extracellular potentials at each node along the fiber. For more info see the documentation on extracellular potentials <extracellular potentials.md> in PyFibers.

loc(loc, target='nodes')

Retrieve a node or section at a given normalized position along the fiber.

Returns either the node or section nearest to loc * (len - 1) from the fiber. Note that the indexing is performed on the nodes or sections list. This means that the number represents the proportion along the list of nodes or sections, not necessarily along the physical fiber (though these are generally the same).

Parameters:
  • loc (float) – Location in the range [0, 1].

  • target (str) – Specifies whether to retrieve from 'nodes' or 'sections'.

Raises:

AssertionError – If loc is not in [0, 1] or if target is not 'nodes' or 'sections'.

Return type:

Section

Returns:

The chosen node or section as a h.Section.

loc_index(loc, target='nodes')

Convert a normalized location [0, 1] into an integer index for nodes or sections.

Parameters:
  • loc (float) – Location in the fiber (from 0 to 1).

  • target (str) – Indicates whether to index into 'nodes' or 'sections'.

Raises:

AssertionError – If loc is not in [0, 1] or if target is invalid.

Return type:

int

Returns:

The integer index corresponding to the node or section.

is_3d()

Check if the fiber is using 3D coordinates.

Return type:

bool

Returns:

True if 3D, otherwise False.

resample_potentials(potentials, potential_coords, center=False, inplace=False)

Use linear interpolation to resample external potentials onto the fiber’s coordinate system (1D).

This is used when extracellular potentials are calculated from an external source, such as a finite element model. The potentials provided by the user should be sampled at high resolution along the fiber’s path and provided alongside the corresponding arc-length coordinates.

If center=True, both the input coordinates and the fiber’s coordinates will be shifted such that their midpoints align.

Parameters:
  • potentials (ndarray) – 1D array of external potential values.

  • potential_coords (ndarray) – 1D array of coordinates corresponding to potentials.

  • center (bool) – If True, center the potentials around the midpoint of each domain.

  • inplace (bool) – If True, update Fiber.potentials with the resampled values.

Return type:

ndarray

Returns:

Interpolated potential values aligned with Fiber.longitudinal_coordinates.

Raises:

AssertionError – If input array sizes or monotonicity checks fail.

apcounts(thresh=-30)

Create NEURON APCount objects at each node to detect action potentials.

Parameters:

thresh (float) – Threshold voltage (mV) for an action potential.

Return type:

None

record_values(ref_attr, allsec=False, indices=None, allow_missing=True, recording_dt=None, recording_tvec=None)

Record NEURON variable references (e.g. membrane voltage) along the fiber.

Note that recording_dt and recording_tvec are mutually exclusive. If both are None, the variable is recorded at every simulation timestep. For more info, see the NEURON docs: https://nrn.readthedocs.io/en/latest/python/programming/math/vector.html#Vector.record

Parameters:
  • ref_attr (str) – The NEURON attribute to record (e.g. '_ref_v').

  • allsec (bool) – If True, record from sections (including nodes). Otherwise, only record from nodes.

  • indices (list[int]) – Specific indices to record from (if None, record from all).

  • allow_missing (bool) – If True, allows missing attributes without raising an error (returns None).

  • recording_dt (float) – The time step [ms] for recording the values (separate from simulation dt). Should be larger than the simulation dt.

  • recording_tvec (h.Vector) –

    A NEURON Vector of time points at which to record the values (ms). Note that the user MUST keep this Vector in memory for the duration of the simulation. This means you must assign it to a variable that is not overwritten or deleted. For example, to record at time points 0, 1, 2, and 3 ms:

    recording_tvec = h.Vector([0, 1, 2, 3])  # store times in a Vector
    fiber.record_values("_ref_v", recording_tvec=recording_tvec)  # pass Vector to record
    stimulation.find_threshold(fiber)  # run the simulation
    plt.plot(recording_tvec, fiber.vm[0])  # plot the recorded values
    

Raises:

ValueError – If indices is an empty list.

Return type:

list[h.Vector | None]

Returns:

A list of NEURON Vector objects or None (if allow_missing=True and the requested attribute is missing).

record_vm(**kwargs)

Record membrane voltage (mV) along the fiber.

Parameters:

kwargs – Additional arguments passed to Fiber.record_values().

Return type:

list[h.Vector | None]

Returns:

A list of NEURON Vector objects recording membrane voltage.

record_im(**kwargs)

Record membrane current (nA) along the fiber.

Parameters:

kwargs – Additional arguments passed to Fiber.record_values().

Return type:

list[h.Vector | None]

Returns:

A list of NEURON Vector objects recording membrane current.

record_vext()

Record extracellular potential (mV) from each section along the fiber.

Return type:

list[Vector]

Returns:

A list of NEURON Vector objects recording extracellular potential.

record_gating(**kwargs)

Record gating parameters (ion channel states) from axon nodes.

The gating variables must be declared in Fiber.gating_variables within the fiber model class.

Parameters:

kwargs – Additional arguments passed to Fiber.record_values().

Return type:

dict[str, list[h.Vector | None]]

Returns:

A dictionary mapping gating variable names to lists of recorded NEURON Vector objects.

Raises:

AssertionError – If Fiber.gating_variables is empty.

set_save_im(**kwargs)

Record membrane current (nA) along the fiber.

Parameters:

kwargs – Additional arguments passed to Fiber.record_values().

Return type:

list[h.Vector | None]

Returns:

A list of NEURON Vector objects recording membrane current.

set_save_vext()

Record extracellular potential (mV) from each section along the fiber.

Return type:

list[Vector]

Returns:

A list of NEURON Vector objects recording extracellular potential.

set_save_vm(**kwargs)

Record membrane voltage (mV) along the fiber.

Parameters:

kwargs – Additional arguments passed to Fiber.record_values().

Return type:

list[h.Vector | None]

Returns:

A list of NEURON Vector objects recording membrane voltage.

set_save_gating(**kwargs)

Record gating parameters (ion channel states) from axon nodes.

The gating variables must be declared in Fiber.gating_variables within the fiber model class.

Parameters:

kwargs – Additional arguments passed to Fiber.record_values().

Return type:

dict[str, list[h.Vector | None]]

Returns:

A dictionary mapping gating variable names to lists of recorded NEURON Vector objects.

Raises:

AssertionError – If Fiber.gating_variables is empty.

point_source_potentials(x, y, z, i0, sigma, inplace=False)

Compute extracellular potentials at the fiber’s coordinates due to a point source stimulus.

See also

Documentation on extracellular potentials <extracellular potentials.md> in PyFibers.

Parameters:
  • x (float) – x-coordinate of the source in µm.

  • y (float) – y-coordinate of the source in µm.

  • z (float) – z-coordinate of the source in µm.

  • i0 (float) – Magnitude of the point-source current (mA).

  • sigma (float | tuple) – Conductivity (S/m). A float for isotropic or a tuple (sigma_x, sigma_y, sigma_z) for anisotropic.

  • inplace (bool) – If True, update Fiber.potentials in-place.

Return type:

ndarray

Returns:

Extracellular potentials at each fiber coordinate, in mV.

measure_cv(start=0.25, end=0.75, tolerance=0.005)

Estimate conduction velocity (m/s) by measuring AP times at two points (start and end).

This method calculates the conduction velocity by comparing the action potential times at two specified normalized locations (using NEURON indexing). It also checks for linear conduction between the two points, within a specified tolerance.

Parameters:
  • start (float) – Starting position for conduction velocity measurement (from 0 to 1).

  • end (float) – Ending position for conduction velocity measurement (from 0 to 1).

  • tolerance (float) – Tolerance (ms) for checking linearity of AP times.

Raises:
  • ValueError – If conduction is not approximately linear between start and end.

  • AssertionError – If no APs are detected at one or both of the measurement nodes.

Return type:

float

Returns:

Conduction velocity in meters per second (m/s).

add_intrinsic_activity(loc=0.1, loc_index=None, avg_interval=1, num_stims=1, start_time=1, noise=0, synapse_tau=0.1, synapse_reversal_potential=0, netcon_weight=0.1)

Add a spike-generating synapse to the fiber for intrinsic activity.

The synapse is generated via a NetStim (spike generator), which is connected to an ExpSyn (exponential synapse) on the chosen node. A NetCon object links them together, injecting an exponentially decaying current upon each spike event.

Parameters:
  • loc (float) – Normalized location along the fiber where the synapse is placed ([0,1]).

  • loc_index (int) – Alternatively, specify an integer index of the node.

  • avg_interval (float) – Average interval between NetStim spikes (ms).

  • num_stims (int) – Number of spikes to deliver.

  • start_time (float) – Time to start delivering spikes (ms).

  • noise (float) – Noise parameter for spike intervals (0 = regular, 1 = Poisson).

  • synapse_tau (float) – Time constant (ms) for synaptic current decay.

  • synapse_reversal_potential (float) – Reversal potential (mV) of the synapse.

  • netcon_weight (float) – Weight of the NetCon between the spike generator and synapse.

Raises:

AssertionError – If neither loc nor loc_index is specified, or if both are specified.

Return type:

None

static calculate_periaxonal_current(from_sec, to_sec, vext_from, vext_to)

Compute the periaxonal current between two compartments for myelinated models.

Parameters:
  • from_sec (Section) – The NEURON h.Section from which current flows.

  • to_sec (Section) – The NEURON h.Section receiving current.

  • vext_from (float) – Extracellular (periaxonal) potential at from_sec (mV).

  • vext_to (float) – Extracellular (periaxonal) potential at to_sec (mV).

Return type:

float

Returns:

The periaxonal current in mA.

membrane_currents(downsample=1)

Compute membrane currents at each section for each time point in the simulation.

Uses the methods described in Pena et. al 2024: http://dx.doi.org/10.1371/journal.pcbi.1011833

For myelinated fibers, this calculation includes periaxonal currents between adjacent sections (based on h.Section.xraxial <neuron.Section.xraxial>). The result is a matrix of shape: [num_timepoints, num_sections].

This method returns a tuple consisting of: 1. A 2D array (time steps x number of sections) of membrane currents in mA. 2. The array of time points corresponding to those currents (downsampled by the specified factor).

Parameters:

downsample (int) – Factor to reduce the temporal resolution (e.g., downsample=2 takes every 2nd time step).

Return type:

ndarray

Returns:

A tuple (i_membrane_matrix, downsampled_time). The matrix contains total currents (mA). The time array contains the corresponding simulation times (ms).

Raises:

RuntimeError – If membrane currents or extracellular potentials were not recorded properly.

static sfap(current_matrix, potentials)

Compute the Single-Fiber Action Potential (SFAP) by multiplying currents with recording potentials.

Parameters:
  • current_matrix (ndarray) – 2D array of shape [timepoints, sections], containing currents in mA.

  • potentials (ndarray) – 1D array of potentials (mV) at each section, length = number of sections.

Raises:

AssertionError – If the number of columns in the current matrix does not match the length of potentials.

Return type:

ndarray

Returns:

The computed SFAP in microvolts (µV).

record_sfap(rec_potentials, downsample=1)

Compute the SFAP time course at a given electrode location.

Parameters:
  • rec_potentials (list | ndarray) – 1D array of precomputed potentials (mV) at each fiber section due to the electrode placement.

  • downsample (int) – Downsampling factor for the time vector (applies to the current matrix).

Return type:

ndarray

Returns:

A tuple (sfap_trace, downsampled_time) where sfap_trace is the computed single-fiber action potential in microvolts (µV) and downsampled_time is the corresponding time array (ms).

set_xyz(x=0, y=0, z=0)

Assign new (x, y, z) shifts to a straight (1D) fiber.

The fiber is assumed to be along the z-axis initially, with x=0 and y=0. This method sets the x and y coordinates to the specified values and shifts the z coordinate by the given amount.

Parameters:
  • x (float) – Shift in the x-direction (µm).

  • y (float) – Shift in the y-direction (µm).

  • z (float) – Shift in the z-direction (µm).

Raises:

AssertionError – If this fiber is a 3D fiber (since this method is for 1D only).

Return type:

None

resample_potentials_3d(potentials, potential_coords, center=False, inplace=False)

Interpolate external potentials onto a 3D fiber coordinate system.

A wrapper around Fiber.resample_potentials() that handles 3D coordinates by first computing the arc length of the provided coordinate array. As with the 1D version, this method is used to resample external potentials (e.g., from a finite element model) onto the fiber.

At present, this does not check that the input coordinates lie exactly along the 3D fiber path. Therefore, it is recommended to use the same coordinates to construct the fiber as to use here. Alternatively, you can create a 1D fiber and calculate the coordinate arc lengths. For more information, see Supplying Extracellular Potentials.

Parameters:
  • potentials (ndarray) – 1D array of external potential values.

  • potential_coords (ndarray) – 2D array of shape (N, 3) representing the (x, y, z) coordinates where the potentials are measured or computed.

  • center (bool) – If True, center the potentials around the midpoint of each domain.

  • inplace (bool) – If True, update Fiber.potentials with the resampled values.

Return type:

ndarray

Returns:

Interpolated potential values aligned with the fiber’s 3D arc-length coordinates.

Raises:

AssertionError – If called on a non-3D fiber or if input coordinate shapes are invalid.

generate(function_list, n_nodes=None, n_sections=None, length=None, enforce_odd_nodecount=True)

Build the fiber model sections in NEURON according to a specified generation strategy.

This method either uses the desired number of nodes (n_nodes), the number of sections (n_sections), or determines the total number of sections from the fiber length and the spacing Fiber.delta_z.

Parameters:
  • function_list (list[Callable]) – List of callable functions that each create a section. Typically structured as [node_func, myelin_func].

  • n_nodes (int) – Number of nodes of Ranvier.

  • n_sections (int) – Total number of sections in the fiber.

  • length (float) – Total length of the fiber (µm). Overrides n_sections if given.

  • enforce_odd_nodecount (bool) – If True, ensures that the fiber has an odd number of nodes.

Return type:

Fiber

Returns:

The updated Fiber instance after generation.

Raises:

AssertionError – If the computed number of sections does not align with the function_list-based pattern.

_create_sections(function_list)

Create and connect NEURON sections for each node or internode in the fiber.

The provided function_list starts with a function for a node, followed by each internodal section in order. Each node is optionally converted to a passive node if it is within the range of Fiber.passive_end_nodes.

Parameters:

function_list (list[Callable]) – A list of functions that each return a new NEURON h.Section.

Return type:

Fiber

Returns:

The updated Fiber instance.

_make_passive(node)

Convert a node section to passive by removing all active mechanisms.

For more info, see Implementations of Fiber Models.

Parameters:

node (Section) – The node h.Section to be made passive.

Return type:

Section

Returns:

The modified section with a passive mechanism inserted.

Raises:

AssertionError – If the node’s name does not contain ‘passive’.

nodebuilder(ind, node_type)

Create a generic node of Ranvier.

This method sets the node length, diameter, and inserts the extracellular mechanism, but does not add active ion channel mechanisms. Subclasses or external code should add the relevant channels to define the node’s electrophysiological behavior.

Parameters:
  • ind (int) – Index of this node in the overall fiber construction.

  • node_type (str) – A string identifying the node as ‘active’ or ‘passive’.

Return type:

Section

Returns:

The newly created node as a NEURON h.Section.