gif_cond_exp_multisynapse#

class brainpy.state.gif_cond_exp_multisynapse(in_size, g_L=Quantity(4., "nS"), E_L=Quantity(-70., "mV"), C_m=Quantity(80., "pF"), V_reset=Quantity(-55., "mV"), Delta_V=Quantity(0.5, "mV"), V_T_star=Quantity(-35., "mV"), lambda_0=1.0, t_ref=Quantity(4., "ms"), tau_syn=(2.0, ), E_rev=(0.0, ), I_e=Quantity(0., "pA"), tau_sfa=(), q_sfa=(), tau_stc=(), q_stc=(), gsl_error_tol=0.001, rng_key=None, V_initializer=Constant(value=-70. mV), spk_fun=ReluGrad(alpha=0.3, width=1.0), spk_reset='hard', name=None)#

Conductance-based generalized integrate-and-fire neuron (GIF) model with multiple synaptic time constants.

gif_cond_exp_multisynapse is the generalized integrate-and-fire neuron according to Mensi et al. (2012) [1] and Pozzorini et al. (2015) [2], with postsynaptic conductances in the form of truncated exponentials and an arbitrary number of synaptic receptor ports.

This is a brainpy.state re-implementation of the NEST simulator model of the same name, using NEST-standard parameterization.

This model features both an adaptation current and a dynamic threshold for spike-frequency adaptation. The membrane potential \(V\) is described by the differential equation:

\[C_\mathrm{m} \frac{dV(t)}{dt} = -g_\mathrm{L}(V(t) - E_\mathrm{L}) - \sum_k g_k(t)(V(t) - E_{\mathrm{rev},k}) - \eta_1(t) - \eta_2(t) - \ldots - \eta_n(t) + I_\mathrm{e} + I_\mathrm{stim}(t)\]

where \(g_k(t)\) are synaptic conductances for receptor port \(k\), each with its own reversal potential \(E_{\mathrm{rev},k}\) and time constant \(\tau_{\mathrm{syn},k}\), and each \(\eta_i\) is a spike-triggered current (stc).

1. Synaptic conductances

Each synaptic conductance decays exponentially:

\[\frac{dg_k}{dt} = -\frac{g_k}{\tau_{\mathrm{syn},k}}\]

On the postsynaptic side, there can be arbitrarily many synaptic time constants. This is achieved by specifying separate receptor ports, each for a different time constant. The number of receptor ports is determined by the length of the tau_syn and E_rev lists, which must have equal length.

2. Spike-triggered currents

Dynamic of each \(\eta_i\) is described by:

\[\tau_{\eta_i} \cdot \frac{d\eta_i}{dt} = -\eta_i\]

and in case of spike emission, its value is increased by a constant:

\[\eta_i = \eta_i + q_{\eta_i} \quad \text{(on spike emission)}\]

3. Spike-frequency adaptation

The neuron produces spikes stochastically according to a point process with the firing intensity:

\[\lambda(t) = \lambda_0 \cdot \exp\left(\frac{V(t) - V_T(t)}{\Delta_V}\right)\]

where \(V_T(t)\) is a time-dependent firing threshold:

\[V_T(t) = V_{T^*} + \gamma_1(t) + \gamma_2(t) + \ldots + \gamma_m(t)\]

where \(\gamma_i\) is a kernel of spike-frequency adaptation (sfa). Dynamic of each \(\gamma_i\) is described by:

\[\tau_{\gamma_i} \cdot \frac{d\gamma_i}{dt} = -\gamma_i\]

and in case of spike emission, its value is increased by a constant:

\[\gamma_i = \gamma_i + q_{\gamma_i} \quad \text{(on spike emission)}\]

4. Stochastic spiking

The probability of firing within a time step \(dt\) is computed using the hazard function:

\[P(\text{spike}) = 1 - \exp(-\lambda(t) \cdot dt)\]

A random number is drawn each (non-refractory) time step and compared to this probability to determine whether a spike occurs.

5. Refractory mechanism

After a spike, the neuron enters an absolute refractory period of duration \(t_\mathrm{ref}\). During this period:

  • \(V_\mathrm{m}\) is clamped to \(V_\mathrm{reset}\),

  • \(dV_\mathrm{m}/dt = 0\),

  • conductances continue to decay,

  • refractory counter decrements each step.

6. Numerical integration and update order

NEST integrates this model with adaptive RKF45. This implementation mirrors that behavior with an RKF45(4,5) integrator and persistent internal step size. The discrete-time update order per simulation step is:

  1. Compute total stc (sum of stc elements) and sfa threshold (V_T_star + sum of sfa elements). Then decay all stc and sfa elements by their respective exponential factors.

  2. Integrate continuous dynamics \([V_\mathrm{m}, g_0, g_1, \ldots, g_{n-1}]\) over \((t, t+dt]\) using RKF45.

  3. Add synaptic conductance jumps from spike inputs arriving this step (per receptor).

  4. If not refractory: compute firing intensity, draw random number, potentially emit spike (update stc/sfa elements, set refractory counter). If refractory: decrement counter, clamp V to V_reset.

  5. Store external current input as \(I_\mathrm{stim}\) for the next step.

7. Multisynapse differences from gif_cond_exp

Unlike gif_cond_exp which has exactly two fixed synaptic channels (excitatory and inhibitory with separate parameters tau_syn_ex, tau_syn_in, E_ex, E_in), this model supports an arbitrary number of receptor ports specified by the lists tau_syn and E_rev.

When connecting to this model, all synaptic weights must be non-negative. Each connection specifies its receptor port via receptor_type, which indexes into the tau_syn / E_rev arrays (1-based in NEST, but 0-based for the add_delta_input interface here).

Note

In the NEST implementation, the stc and sfa element jumps occur immediately after spike emission. The GIF toolbox uses a different convention where jumps occur after the refractory period. Conversion:

\[q_{\eta,\text{toolbox}} = q_{\eta,\text{NEST}} \cdot (1 - \exp(-t_\mathrm{ref} / \tau_\eta))\]

Note

Because spiking is stochastic (random number drawn each step), exact spike-time reproducibility requires matching the random number generator state. For deterministic testing, set rng_key explicitly.

Parameters:
  • in_size (int, sequence of int) – Population shape (e.g., 100 or (10, 10)). Required.

  • g_L (ArrayLike, default: 4.0 nS) – Leak conductance. Must be strictly positive. Shape: scalar or broadcastable to in_size.

  • E_L (ArrayLike, default: -70.0 mV) – Leak reversal potential (resting potential). Shape: scalar or broadcastable to in_size.

  • C_m (ArrayLike, default: 80.0 pF) – Membrane capacitance. Must be strictly positive. Shape: scalar or broadcastable to in_size.

  • V_reset (ArrayLike, default: -55.0 mV) – Reset potential after spike. Shape: scalar or broadcastable to in_size.

  • Delta_V (ArrayLike, default: 0.5 mV) – Stochasticity level for exponential firing intensity. Must be strictly positive. Shape: scalar or broadcastable to in_size.

  • V_T_star (ArrayLike, default: -35.0 mV) – Base (non-adapting) firing threshold. Shape: scalar or broadcastable to in_size.

  • lambda_0 (float, default: 1.0) – Stochastic intensity at threshold (in 1/s). Must be non-negative. Internally converted to 1/ms.

  • t_ref (ArrayLike, default: 4.0 ms) – Absolute refractory period duration. Must be non-negative. Shape: scalar or broadcastable to in_size.

  • tau_syn (Sequence[float], default: (2.0,)) – Synaptic conductance time constants for each receptor port (in ms). Each element must be strictly positive. Must have same length as E_rev. At least one element required.

  • E_rev (Sequence[float], default: (0.0,)) – Reversal potentials for each receptor port (in mV). Must have same length as tau_syn. At least one element required.

  • I_e (ArrayLike, default: 0.0 pA) – Constant external current. Shape: scalar or broadcastable to in_size.

  • tau_sfa (Sequence[float], default: ()) – Time constants for spike-frequency adaptation (SFA) threshold elements (in ms). Each element must be strictly positive. Must have same length as q_sfa.

  • q_sfa (Sequence[float], default: ()) – Jump values for SFA threshold elements (in mV). Must have same length as tau_sfa.

  • tau_stc (Sequence[float], default: ()) – Time constants for spike-triggered current (STC) elements (in ms). Each element must be strictly positive. Must have same length as q_stc.

  • q_stc (Sequence[float], default: ()) – Jump values for STC elements (in nA). Must have same length as tau_stc.

  • gsl_error_tol (ArrayLike, default: 1e-3) – Unitless local RKF45 error tolerance, broadcastable and strictly positive.

  • rng_key (jax.Array, optional) – JAX PRNG key for stochastic spiking. If None, defaults to jax.random.PRNGKey(0).

  • V_initializer (Callable, default: Constant(-70.0 mV)) – Initializer for membrane potential. Must return values compatible with in_size.

  • spk_fun (Callable, default: ReluGrad()) – Surrogate gradient function for spike generation. Used in gradient-based learning.

  • spk_reset (str, default: 'hard') – Spike reset mode. ‘hard’ (stop gradient, matches NEST) or ‘soft’ (subtract threshold).

  • name (str, optional) – Module name. If None, auto-generated.

Parameter Mapping

Maps brainpy.state parameter names to NEST equivalents for cross-framework compatibility:

Parameter

Default

Math equivalent

Description

in_size

(required)

Population shape

g_L

4.0 nS

\(g_\mathrm{L}\)

Leak conductance

E_L

-70.0 mV

\(E_\mathrm{L}\)

Leak reversal potential

C_m

80.0 pF

\(C_\mathrm{m}\)

Membrane capacitance

V_reset

-55.0 mV

\(V_\mathrm{reset}\)

Reset potential

Delta_V

0.5 mV

\(\Delta_V\)

Stochasticity level

V_T_star

-35.0 mV

\(V_{T^*}\)

Base firing threshold

lambda_0

1.0 /s

\(\lambda_0\)

Stochastic intensity at threshold

t_ref

4.0 ms

\(t_\mathrm{ref}\)

Absolute refractory period

tau_syn

(2.0,) ms

\(\tau_{\mathrm{syn},k}\)

Synaptic conductance time constants (list)

E_rev

(0.0,) mV

\(E_{\mathrm{rev},k}\)

Reversal potentials (list, same size as tau_syn)

I_e

0.0 pA

\(I_\mathrm{e}\)

Constant external current

tau_sfa

() ms

\(\tau_{\gamma_i}\)

SFA time constants (tuple/list)

q_sfa

() mV

\(q_{\gamma_i}\)

SFA jump values (tuple/list)

tau_stc

() ms

\(\tau_{\eta_i}\)

STC time constants (tuple/list)

q_stc

() nA

\(q_{\eta_i}\)

STC jump values (tuple/list)

gsl_error_tol

1e-3

Local absolute tolerance for RKF45 error estimate

rng_key

None

JAX PRNG key for stochastic spiking

V_initializer

Constant(-70 mV)

Initializer for membrane potential

spk_fun

ReluGrad()

Surrogate spike function

spk_reset

'hard'

Reset mode; hard reset matches NEST

State Variables

After init_state(), the following state variables are available:

State variable

Type

Description

V

HiddenState

Membrane potential \(V_\mathrm{m}\) (mV)

g

List[HiddenState]

List of synaptic conductances \(g_k\) (nS), one per receptor

refractory_step_count

ShortTermState

Remaining refractory grid steps (int32)

integration_step

ShortTermState

Internal RKF45 step-size state (ms)

I_stim

ShortTermState

Buffered current applied in next step (pA)

last_spike_time

ShortTermState

Last spike time (ms)

Additionally, the following NumPy arrays are maintained internally:

  • _stc_elems – shape (len(tau_stc), *in_size) – individual stc elements (nA)

  • _sfa_elems – shape (len(tau_sfa), *in_size) – individual sfa elements (mV)

  • _stc_val – shape in_size – total spike-triggered current (nA)

  • _sfa_val – shape in_size – adaptive threshold \(V_T(t)\) (mV)

Raises:

ValueError – If C_m <= 0, g_L <= 0, Delta_V <= 0, t_ref < 0, lambda_0 < 0, any tau_syn <= 0, any tau_sfa <= 0, any tau_stc <= 0, len(tau_syn) != len(E_rev), len(tau_syn) == 0, len(tau_sfa) != len(q_sfa), or len(tau_stc) != len(q_stc).

Notes

  • Defaults follow NEST C++ source for gif_cond_exp_multisynapse.

  • lambda_0 is specified in 1/s (as in NEST’s Python interface) and is internally converted to 1/ms for computation.

  • All synaptic spike weights must be non-negative (conductance-based model).

  • Delta inputs for synaptic conductances are indexed by receptor port (0-based). Use add_delta_input(f'receptor_{port}', weight * u.nS) to add a conductance jump to a specific receptor port. If no receptor port is specified in the key, the input defaults to receptor 0.

  • RKF45 integration with adaptive step size ensures numerical stability for stiff systems, matching NEST’s GSL-based integrator behavior.

  • The stochastic spiking mechanism uses JAX PRNG, which is split each time step to ensure reproducible randomness under JIT compilation.

  • The model supports an arbitrary number of receptor ports, making it suitable for modeling neurons with multiple synaptic receptor types (e.g., AMPA, NMDA, GABA_A, GABA_B).

References

See also

gif_cond_exp

Two-channel GIF model with fixed excitatory/inhibitory receptors

iaf_cond_exp

Conductance-based leaky integrate-and-fire

gif_psc_exp_multisynapse

Current-based GIF with multiple synaptic time constants

Examples

Create a GIF neuron with three receptor ports (e.g., AMPA, NMDA, GABA_A):

>>> import brainpy.state as bst
>>> import saiunit as u
>>> # AMPA: fast excitatory, NMDA: slow excitatory, GABA_A: fast inhibitory
>>> neuron = bst.gif_cond_exp_multisynapse(
...     in_size=100,
...     tau_syn=(2.0, 20.0, 5.0),  # ms
...     E_rev=(0.0, 0.0, -85.0),   # mV
...     tau_sfa=(100.0,),           # ms
...     q_sfa=(5.0,),               # mV
...     tau_stc=(50.0,),            # ms
...     q_stc=(10.0,)               # nA
... )
>>> neuron.init_all_states()
>>> # Add AMPA input to receptor 0
>>> neuron.add_delta_input('receptor_0', 0.5 * u.nS)
>>> # Add NMDA input to receptor 1
>>> neuron.add_delta_input('receptor_1', 0.3 * u.nS)
>>> # Add GABA_A input to receptor 2
>>> neuron.add_delta_input('receptor_2', 0.8 * u.nS)
get_spike(V=None)[source]#

Compute spike output using surrogate gradient function.

This method is used for gradient-based learning. The actual spike emission during simulation is stochastic and handled in update().

Parameters:

V (ArrayLike, optional) – Membrane potential. If None, uses self.V.value. Shape: (*batch_size, *in_size) with unit mV.

Returns:

Spike values in [0, 1] range. Shape matches input V. Dimensionless. Uses spk_fun to compute differentiable spike output from scaled voltage.

Return type:

ArrayLike

Notes

The spike output is computed as spk_fun((V - V_reset) / Delta_V), providing a differentiable approximation for gradient-based learning. This is separate from the stochastic spike mechanism used in update().

init_state(**kwargs)[source]#

Initialize all state variables.

Parameters:

**kwargs – Unused compatibility parameters accepted by the base-state API.

Notes

This method initializes: - V: membrane potential to V_initializer values - g: list of n_receptors conductance states, all initialized to 0 nS - last_spike_time: initialized to -1e7 ms (far past) - refractory_step_count: initialized to 0 - integration_step: initialized to current dt - I_stim: initialized to 0 pA - _stc_elems, _sfa_elems: NumPy arrays for adaptation elements, all zeros - _stc_val, _sfa_val: computed totals - _rng_state: JAX PRNG state

property n_receptors#

Number of synaptic receptor ports.

Returns:

Number of receptor ports, equal to len(tau_syn) and len(E_rev).

Return type:

int

reset_state(batch_size=None, **kwargs)[source]#

Reset all state variables to their initial values.

Resets membrane potential, conductances, refractory counter, integration step, stimulus buffer, adaptation elements, and RNG state.

Parameters:
  • batch_size (int, optional) – Unused; present for API compatibility.

  • **kwargs – Unused compatibility parameters.

update(x=Quantity(0., 'pA'))[source]#

Update neuron state for one time step.

Performs the following operations in order: 1. Decay adaptation elements (stc, sfa) and compute totals 2. Integrate membrane and conductance dynamics using RKF45 3. Add synaptic conductance jumps from spike inputs 4. Stochastic spike check (if not refractory) or refractory countdown 5. Store external current for next step

Parameters:

x (ArrayLike, default: 0.0 pA) – External current input for this time step. Shape: scalar or broadcastable to in_size. Unit: pA. This is stored as I_stim for use in the next time step.

Returns:

Binary spike output (0 or 1) for this time step. Shape: (*batch_size, *in_size). dtype: float32. Value 1 indicates spike emission.

Return type:

jnp.ndarray

Notes

The update implements NEST’s exact algorithm:

  • Adaptation elements decay exponentially by factor exp(-dt/tau)

  • RKF45 integration uses adaptive step size (stored in integration_step state)

  • Spike probability computed as 1 - exp(-lambda * dt) where lambda = lambda_0 * exp((V - V_T) / Delta_V)

  • Random number generation uses JAX PRNG with automatic state splitting

  • Refractory neurons have V clamped to V_reset and cannot spike

  • All computations performed element-wise in NumPy for each neuron

The current input x is buffered (not used in current step) to match NEST’s one-step delay convention.