aeif_cond_beta_multisynapse#

class brainpy.state.aeif_cond_beta_multisynapse(in_size, V_peak=Quantity(0., "mV"), V_reset=Quantity(-60., "mV"), t_ref=Quantity(0., "ms"), g_L=Quantity(30., "nS"), C_m=Quantity(281., "pF"), E_L=Quantity(-70.6, "mV"), Delta_T=Quantity(2., "mV"), tau_w=Quantity(144., "ms"), a=Quantity(4., "nS"), b=Quantity(80.5, "pA"), V_th=Quantity(-50.4, "mV"), tau_rise=Quantity([2.], "ms"), tau_decay=Quantity([20.], "ms"), E_rev=Quantity([0.], "mV"), I_e=Quantity(0., "pA"), gsl_error_tol=1e-06, V_initializer=Constant(value=-70.6 mV), g_initializer=Constant(value=0. nS), w_initializer=Constant(value=0. pA), spk_fun=ReluGrad(alpha=0.3, width=1.0), spk_reset='hard', ref_var=False, name=None)#

NEST-compatible aeif_cond_beta_multisynapse neuron model.

Conductance-based adaptive exponential integrate-and-fire neuron with beta-shaped synaptic conductances and an arbitrary number of receptor ports. Implements NEST’s aeif_cond_beta_multisynapse model with source-level parity in update ordering, refractory handling, and spike detection.

This model extends the adaptive exponential integrate-and-fire (AdEx) framework [1] with beta-function synaptic conductances instead of exponential or alpha shapes. Each receptor port maintains independent rise/decay time constants and reversal potentials, enabling multi-receptor networks (e.g., AMPA + GABA_A).

Parameters:
  • in_size (Size) – Population shape as integer, tuple, or Size object. Required.

  • V_peak (ArrayLike, optional) – Spike detection threshold (mV). Used when Delta_T > 0; otherwise V_th is used. Must satisfy V_peak >= V_th. Default: 0.0 mV.

  • V_reset (ArrayLike, optional) – Post-spike reset potential (mV). Must satisfy V_reset < V_peak. Default: -60.0 mV.

  • t_ref (ArrayLike, optional) – Absolute refractory period (ms). During refractory, dV/dt = 0 and voltage is clamped to V_reset. Default: 0.0 ms (no refractory).

  • g_L (ArrayLike, optional) – Leak conductance (nS). Must be strictly positive. Default: 30.0 nS.

  • C_m (ArrayLike, optional) – Membrane capacitance (pF). Must be strictly positive. Default: 281.0 pF.

  • E_L (ArrayLike, optional) – Leak reversal potential (mV). Default: -70.6 mV.

  • Delta_T (ArrayLike, optional) – Exponential slope factor (mV). Must be non-negative. When Delta_T = 0, reduces to LIF-like dynamics. Default: 2.0 mV.

  • tau_w (ArrayLike, optional) – Adaptation time constant (ms). Must be strictly positive. Default: 144.0 ms.

  • a (ArrayLike, optional) – Subthreshold adaptation coupling (nS). Default: 4.0 nS.

  • b (ArrayLike, optional) – Spike-triggered adaptation increment (pA). Added to w on each spike. Default: 80.5 pA.

  • V_th (ArrayLike, optional) – Spike initiation threshold (mV) for exponential term. Must satisfy V_th <= V_peak. Default: -50.4 mV.

  • tau_rise (ArrayLike, optional) – Synaptic rise time constants (ms) per receptor, shape (n_receptors,). Must be strictly positive and satisfy tau_rise <= tau_decay element-wise. Default: (2.0,) ms (single receptor).

  • tau_decay (ArrayLike, optional) – Synaptic decay time constants (ms) per receptor, shape (n_receptors,). Must be strictly positive and satisfy tau_decay >= tau_rise element-wise. Default: (20.0,) ms (single receptor).

  • E_rev (ArrayLike, optional) – Reversal potentials (mV) per receptor, shape (n_receptors,). Default: (0.0,) mV (excitatory-like).

  • I_e (ArrayLike, optional) – Constant external current (pA). Default: 0.0 pA.

  • gsl_error_tol (ArrayLike, optional) – RKF45 local error tolerance (unitless). Smaller values improve accuracy but increase computational cost. Must be strictly positive. Default: 1e-6.

  • V_initializer (Callable, optional) – Membrane potential initializer. Default: Constant(-70.6 mV).

  • g_initializer (Callable, optional) – Conductance state initializer with shape [..., n_receptors]. Default: Constant(0.0 nS).

  • w_initializer (Callable, optional) – Adaptation current initializer. Default: Constant(0.0 pA).

  • spk_fun (Callable, optional) – Surrogate gradient function for differentiable spike generation. Default: ReluGrad().

  • spk_reset (str, optional) – Spike reset mode. 'hard' (stop_gradient, matches NEST) or 'soft' (subtract threshold). Default: 'hard'.

  • ref_var (bool, optional) – If True, expose refractory state variable (boolean indicator). Default: False.

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

Parameter Mapping

BrainPy Parameter

NEST Parameter

Description

in_size

(model count)

Population shape

V_peak

V_peak

Spike detection threshold

V_reset

V_reset

Reset potential

t_ref

t_ref

Refractory period

g_L

g_L

Leak conductance

C_m

C_m

Membrane capacitance

E_L

E_L

Leak reversal

Delta_T

Delta_T

Slope factor

tau_w

tau_w

Adaptation time constant

a

a

Subthreshold adaptation

b

b

Spike-triggered adaptation

V_th

V_th

Exponential threshold

tau_rise

tau_rise

Rise time per receptor

tau_decay

tau_decay

Decay time per receptor

E_rev

E_rev

Reversal potential per receptor

I_e

I_e

Constant external current

gsl_error_tol

gsl_error_tol

RKF45 tolerance

Mathematical Model

1. Membrane Dynamics

The membrane voltage \(V\) evolves according to:

\[C_m \frac{dV}{dt} = -g_L (V - E_L) + g_L \Delta_T \exp\!\left(\frac{V - V_{th}}{\Delta_T}\right) + \sum_{k=1}^{n_{\text{rec}}} g_k (E_{\text{rev},k} - V) - w + I_e + I_{\text{stim}},\]

where:

  • \(C_m\) – membrane capacitance (pF)

  • \(g_L\) – leak conductance (nS)

  • \(E_L\) – leak reversal potential (mV)

  • \(\Delta_T\) – exponential slope factor (mV)

  • \(V_{th}\) – spike initiation threshold (mV)

  • \(g_k\) – conductance of receptor \(k\) (nS)

  • \(E_{\text{rev},k}\) – reversal potential of receptor \(k\) (mV)

  • \(w\) – adaptation current (pA)

  • \(I_e\) – constant external current (pA)

  • \(I_{\text{stim}}\) – delayed injected current (pA)

During refractory period, \(dV/dt = 0\) and \(V\) is clamped to \(V_{\text{reset}}\). Outside refractory, the exponential term uses \(\min(V, V_{\text{peak}})\) to prevent numerical overflow.

2. Adaptation Dynamics

The adaptation current \(w\) follows:

\[\tau_w \frac{dw}{dt} = a (V - E_L) - w,\]

where \(a\) (nS) couples subthreshold membrane voltage fluctuations to adaptation. On each spike, \(w \leftarrow w + b\) implements spike- triggered adaptation.

3. Beta-Function Synaptic Conductances

Each receptor \(k\) maintains two state variables:

\[\frac{d\,dg_k}{dt} = -\frac{dg_k}{\tau_{\text{rise},k}}, \quad \frac{dg_k}{dt} = dg_k - \frac{g_k}{\tau_{\text{decay},k}}.\]

Incoming spikes with weight \(w_k\) (nS) increment the auxiliary state:

\[dg_k \leftarrow dg_k + g_{0,k} w_k,\]

where \(g_{0,k}\) is the beta normalization factor ensuring unit weight produces unit peak conductance:

\[g_{0,k} = \frac{1/\tau_{\text{rise},k} - 1/\tau_{\text{decay},k}}{\exp(-t_{\text{peak}}/\tau_{\text{decay},k}) - \exp(-t_{\text{peak}}/\tau_{\text{rise},k})},\]

with \(t_{\text{peak}} = \tau_{\text{decay},k} \tau_{\text{rise},k} \log(\tau_{\text{decay},k}/\tau_{\text{rise},k}) / (\tau_{\text{decay},k} - \tau_{\text{rise},k})\). In the equal-time-constant limit, this reduces to the alpha normalization \(e / \tau\).

4. Spike Detection and Reset

A spike is detected when:

  • \(V \geq V_{\text{peak}}\) if \(\Delta_T > 0\)

  • \(V \geq V_{th}\) if \(\Delta_T = 0\)

Upon spike detection (within RKF45 substeps):

  1. \(V \leftarrow V_{\text{reset}}\)

  2. \(w \leftarrow w + b\)

  3. Refractory counter \(r \leftarrow \lceil t_{\text{ref}} / dt \rceil + 1\) (if t_ref > 0)

5. Update Order (NEST Semantics)

Each simulation step \((t, t+dt]\) proceeds as:

  1. Integrate ODEs using adaptive RKF45 with internal substeps

  2. Inside integration: apply refractory clamp and spike/reset logic

  3. Decrement refractory counter once (outside integration)

  4. Apply incoming spike events to dg states with beta normalization

  5. Store continuous current input for next step (one-step delay)

Computational Notes

  • Numerical integration: Runge-Kutta-Fehlberg (RKF45) adaptive solver with local error tolerance gsl_error_tol. Internal step size adapts dynamically and persists across simulation steps.

  • Refractory handling: During refractory, effective voltage is clamped to V_reset for all ODE terms, including adaptation coupling.

  • Overflow protection: Exponential term uses \(\min(V, V_{\text{peak}})\) outside refractory to prevent \(\exp(\cdot)\) overflow. Validation ensures \((V_{\text{peak}} - V_{th}) / \Delta_T\) stays below overflow threshold when \(\Delta_T > 0\).

  • Spike event format: spike_events must be an iterable of (receptor_type, weight) tuples or dicts with keys receptor_type/ receptor and weight. Receptor types are 1-based (NEST convention): 1 <= receptor_type <= n_receptors. Weights (nS) must be non-negative.

  • Default input mapping: add_delta_input stream is mapped to receptor 1; weights must be non-negative.

  • Instability detection: Integration raises ValueError if \(V < -1000\) mV or \(|w| > 10^6\) pA, indicating numerical collapse.

V#

Membrane potential (mV), shape (*in_size,).

Type:

HiddenState

w#

Adaptation current (pA), shape (*in_size,).

Type:

HiddenState

dg#

Beta auxiliary states (nS/ms), shape (*in_size, n_receptors).

Type:

ShortTermState

g#

Receptor conductances (nS), shape (*in_size, n_receptors).

Type:

HiddenState

refractory_step_count#

Remaining refractory steps (int32), shape (*in_size,).

Type:

ShortTermState

integration_step#

Persistent RKF45 step size (ms), shape (*in_size,).

Type:

ShortTermState

I_stim#

One-step delayed current buffer (pA), shape (*in_size,).

Type:

ShortTermState

last_spike_time#

Last spike time (ms), shape (*in_size,). Initialized to -1e7 ms.

Type:

ShortTermState

refractory#

Boolean refractory indicator, shape (*in_size,). Only present if ref_var=True.

Type:

ShortTermState, optional

Raises:
  • ValueError – If tau_rise.size != tau_decay.size != E_rev.size.

  • ValueError – If any tau_rise <= 0 or tau_decay <= 0.

  • ValueError – If any tau_decay < tau_rise.

  • ValueError – If any V_peak < V_th or V_reset >= V_peak.

  • ValueError – If Delta_T < 0 or C_m <= 0 or t_ref < 0 or tau_w <= 0.

  • ValueError – If gsl_error_tol <= 0.

  • ValueError – If \((V_{\text{peak}} - V_{th}) / \Delta_T\) exceeds overflow threshold (when \(\Delta_T > 0\)).

  • ValueError – During update, if receptor type out of range [1, n_receptors].

  • ValueError – During update, if synaptic weight is negative (conductance constraint).

  • ValueError – During update, if numerical instability detected (\(V < -1000\) mV or \(|w| > 10^6\) pA).

See also

aeif_cond_alpha_multisynapse

Alpha-function variant

aeif_cond_exp

Single exponential synapse

aeif_psc_exp

Current-based AdEx

Notes

  • Default t_ref = 0 matches NEST and allows multiple spikes per timestep. Set t_ref > 0 to enforce physiological refractory periods.

  • Beta conductances provide more realistic synaptic shapes than single exponentials but require two state variables per receptor (dg and g).

  • When tau_rise = tau_decay, normalization degenerates to alpha-function limit \(e / \tau\).

References

Examples

Create a two-receptor neuron (excitatory + inhibitory):

>>> import brainpy.state as bp
>>> import saiunit as u
>>> neuron = bp.aeif_cond_beta_multisynapse(
...     in_size=10,
...     tau_rise=(2.0, 0.5) * u.ms,
...     tau_decay=(20.0, 8.0) * u.ms,
...     E_rev=(0.0, -80.0) * u.mV,  # excitatory, inhibitory
... )
>>> neuron.n_receptors
2

Simulate with receptor-specific spike events:

>>> import brainstate as bst
>>> with bst.environ.context(dt=0.1 * u.ms):
...     neuron.init_all_states()
...     # Excitatory spike to receptor 1
...     events = [(1, 5.0 * u.nS)]
...     spk = neuron.update(x=100.0 * u.pA, spike_events=events)
...     print(neuron.V.value)

Multi-receptor dictionary format:

>>> events = [
...     {'receptor_type': 1, 'weight': 3.0 * u.nS},
...     {'receptor_type': 2, 'weight': 2.0 * u.nS},
... ]
>>> spk = neuron.update(spike_events=events)
__init__(in_size, V_peak=Quantity(0., "mV"), V_reset=Quantity(-60., "mV"), t_ref=Quantity(0., "ms"), g_L=Quantity(30., "nS"), C_m=Quantity(281., "pF"), E_L=Quantity(-70.6, "mV"), Delta_T=Quantity(2., "mV"), tau_w=Quantity(144., "ms"), a=Quantity(4., "nS"), b=Quantity(80.5, "pA"), V_th=Quantity(-50.4, "mV"), tau_rise=Quantity([2.], "ms"), tau_decay=Quantity([20.], "ms"), E_rev=Quantity([0.], "mV"), I_e=Quantity(0., "pA"), gsl_error_tol=1e-06, V_initializer=Constant(value=-70.6 mV), g_initializer=Constant(value=0. nS), w_initializer=Constant(value=0. pA), spk_fun=ReluGrad(alpha=0.3, width=1.0), spk_reset='hard', ref_var=False, name=None)[source]#

Initialize aeif_cond_beta_multisynapse neuron.

All parameters are documented in the class docstring.

get_spike(V=None)[source]#

Compute surrogate spike output for gradient-based learning.

Applies surrogate gradient function to scaled membrane potential for differentiable spike generation. Does not modify state variables.

Parameters:

V (ArrayLike, optional) – Membrane potential (mV). If None, uses current self.V.value.

Returns:

Surrogate spike output in [0, 1], shape (*in_size,). Produced by spk_fun applied to (V - V_th) / (V_th - V_reset).

Return type:

ArrayLike

Notes

This method is primarily used for gradient computation in training contexts. Actual spike detection during forward simulation uses hard thresholds in update method.

init_state(**kwargs)[source]#

Initialize all state variables.

Creates HiddenState and ShortTermState attributes for membrane potential, adaptation current, receptor conductances, refractory counters, RKF45 step size, and delayed current buffer.

Parameters:

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

Notes

Initializes: - V (HiddenState): membrane potential from V_initializer - w (HiddenState): adaptation current from w_initializer - dg (ShortTermState): beta auxiliary states, initialized to zero - g (HiddenState): receptor conductances from g_initializer - last_spike_time (ShortTermState): initialized to -1e7 ms - refractory_step_count (ShortTermState): initialized to 0 - integration_step (ShortTermState): RKF45 step size, initialized to dt - I_stim (ShortTermState): delayed current buffer, initialized to 0 pA - refractory (ShortTermState, optional): boolean indicator if ref_var=True

property n_receptors#

Number of receptor ports.

Returns:

Number of receptor types, inferred from tau_rise.size.

Return type:

int

property recordables#

List of recordable state variable names.

Returns:

Dynamic recordables following NEST naming: ['V_m', 'w', 'g_1', 'g_2', ..., 'g_n'].

Return type:

list of str

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

Advance model by one simulation timestep (NEST-compatible update).

Integrates ODEs over \((t, t+dt]\) using adaptive RKF45 with vectorized integration, spike detection, refractory handling, and receptor-specific spike event application. Follows NEST’s update ordering exactly.

Parameters:
  • x (ArrayLike, optional) – Continuous current input (pA), shape broadcastable to (*in_size,). Summed with current_inputs and I_e, then delayed by one timestep (NEST semantics). Default: 0.0 pA.

  • spike_events (Iterable or None, optional) – Incoming spike events as: - List of (receptor_type, weight) tuples - List of dicts with keys 'receptor_type'/'receptor' and 'weight' - Single dict (auto-wrapped to list) - None (no spike input) Receptor types are 1-based: 1 <= receptor_type <= n_receptors. Weights (nS) must be non-negative. Default: None.

Returns:

Binary spike indicator (0 or 1), shape (*in_size,). Float64 for gradient compatibility. Value is 1.0 if spike occurred during \((t, t+dt]\), else 0.0.

Return type:

ArrayLike

Raises:
  • ValueError – If receptor type out of range [1, n_receptors].

  • ValueError – If any spike event weight is negative (conductance constraint).

  • ValueError – If add_delta_input stream contains negative values (mapped to receptor 1).

  • ValueError – If no receptor ports exist but delta_inputs or spike_events are non-zero.

  • ValueError – If numerical instability detected during integration (\(V < -1000\) mV or \(|w| > 10^6\) pA).

Notes

Integration is performed with an adaptive vectorized RKF45 loop, including in-loop spike/reset/adaptation events and optional multiple spikes per step. All arithmetic is unit-aware via saiunit.math.