aeif_psc_exp#
- class brainpy.state.aeif_psc_exp(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_syn_ex=Quantity(0.2, "ms"), tau_syn_in=Quantity(2., "ms"), I_e=Quantity(0., "pA"), gsl_error_tol=1e-06, V_initializer=Constant(value=-70.6 mV), I_ex_initializer=Constant(value=0. pA), I_in_initializer=Constant(value=0. pA), 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 adaptive exponential integrate-and-fire neuron with exponential synapses.
Current-based adaptive exponential integrate-and-fire neuron with exponentially decaying synaptic currents. Implements the AdEx model of Brette & Gerstner (2005) with spike-triggered adaptation, subthreshold adaptation coupling, and separate excitatory/inhibitory exponential current synapses. Follows NEST
models/aeif_psc_exp.{h,cpp}implementation exactly.1. Mathematical Model
Membrane and adaptation dynamics:
The membrane potential \(V\) and adaptation current \(w\) evolve as:
\[C_m \frac{dV}{dt} = -g_L (V - E_L) + g_L \Delta_T \exp\!\left(\frac{V - V_{th}}{\Delta_T}\right) - w + I_{ex} - I_{in} + I_e + I_{stim}\]\[\tau_w \frac{dw}{dt} = a (V - E_L) - w\]where \(C_m\) is membrane capacitance, \(g_L\) is leak conductance, \(E_L\) is leak reversal, \(\Delta_T\) is the exponential slope factor, \(V_{th}\) is the spike threshold, \(a\) couples subthreshold voltage to adaptation, and \(\tau_w\) is the adaptation time constant.
Synaptic current dynamics:
Excitatory and inhibitory currents decay exponentially:
\[\frac{d I_{ex}}{dt} = -\frac{I_{ex}}{\tau_{syn,ex}}, \qquad \frac{d I_{in}}{dt} = -\frac{I_{in}}{\tau_{syn,in}}\]Incoming spike weights (in pA) are split by sign and applied instantaneously:
\[I_{ex} \leftarrow I_{ex} + \max(w, 0), \qquad I_{in} \leftarrow I_{in} + \max(-w, 0)\]2. Refractory and Spike Handling (NEST Semantics)
During refractory period (\(r > 0\) steps remaining), the effective voltage used in the RHS is clamped to \(V_{\text{reset}}\) and \(dV/dt = 0\). Outside refractory, the effective voltage is \(\min(V, V_{\text{peak}})\).
- Spike detection threshold:
\(V_{\text{peak}}\) if \(\Delta_T > 0\) (exponential regime)
\(V_{th}\) if \(\Delta_T = 0\) (integrate-and-fire limit)
- On each detected spike:
\(V \leftarrow V_{\text{reset}}\)
\(w \leftarrow w + b\) (spike-triggered adaptation increment)
Refractory counter set to
refractory_counts + 1(ift_ref > 0)
Spike detection/reset occurs inside the RKF45 substep loop. With
t_ref = 0, multiple spikes can occur within one simulation step, matching NEST behavior.3. Update Order Per Simulation Step
Integrate ODEs on \((t, t+dt]\) via adaptive RKF45 (Runge-Kutta-Fehlberg 4(5))
Inside integration loop: apply refractory clamp, detect spike, reset, adapt
After integration: decrement refractory counter by 1
Apply arriving spike weights to \(I_{ex}\), \(I_{in}\)
Store external current input \(x\) into one-step delayed buffer \(I_{\text{stim}}\)
4. Numerical Integration
Uses adaptive RKF45 with local error control. Step size \(h\) is adjusted to keep error below
gsl_error_tol. Integration step size is persistent across simulation steps for efficiency.- Parameters:
in_size (
intortupleofint) – Population shape. Scalar for 1D population, tuple for multi-dimensional.V_peak (
ArrayLike, optional) – Spike detection threshold (ifDelta_T > 0). Units: mV. Default: 0.0 mV. Scalar or broadcastable toin_size.V_reset (
ArrayLike, optional) – Reset potential after spike. Units: mV. Default: -60.0 mV. Scalar or broadcastable toin_size. Must satisfyV_reset < V_peak.t_ref (
ArrayLike, optional) – Absolute refractory period duration. Units: ms. Default: 0.0 ms. Scalar or broadcastable toin_size. Zero allows multiple spikes per step.g_L (
ArrayLike, optional) – Leak conductance. Units: nS. Default: 30.0 nS. Scalar or broadcastable toin_size. Must be positive.C_m (
ArrayLike, optional) – Membrane capacitance. Units: pF. Default: 281.0 pF. Scalar or broadcastable toin_size. Must be positive.E_L (
ArrayLike, optional) – Leak reversal potential. Units: mV. Default: -70.6 mV. Scalar or broadcastable toin_size.Delta_T (
ArrayLike, optional) – Exponential slope factor. Units: mV. Default: 2.0 mV. Scalar or broadcastable toin_size. Zero recovers integrate-and-fire limit. Must be non-negative. Large values relative toV_peak - V_thmay cause overflow.tau_w (
ArrayLike, optional) – Adaptation time constant. Units: ms. Default: 144.0 ms. Scalar or broadcastable toin_size. Must be positive.a (
ArrayLike, optional) – Subthreshold adaptation coupling. Units: nS. Default: 4.0 nS. Scalar or broadcastable toin_size. Couples voltage deviation to adaptation.b (
ArrayLike, optional) – Spike-triggered adaptation increment. Units: pA. Default: 80.5 pA. Scalar or broadcastable toin_size. Added towon each spike.V_th (
ArrayLike, optional) – Spike initiation threshold (in exponential term). Units: mV. Default: -50.4 mV. Scalar or broadcastable toin_size. Must satisfyV_th <= V_peak.tau_syn_ex (
ArrayLike, optional) – Excitatory synaptic current time constant. Units: ms. Default: 0.2 ms. Scalar or broadcastable toin_size. Must be positive.tau_syn_in (
ArrayLike, optional) – Inhibitory synaptic current time constant. Units: ms. Default: 2.0 ms. Scalar or broadcastable toin_size. Must be positive.I_e (
ArrayLike, optional) – Constant external current. Units: pA. Default: 0.0 pA. Scalar or broadcastable toin_size.gsl_error_tol (
ArrayLike, optional) – RKF45 local error tolerance. Dimensionless. Default: 1e-6. Scalar or broadcastable toin_size. Must be positive. Smaller values increase accuracy at the cost of smaller integration steps.V_initializer (
Callable, optional) – Membrane potential initializer. Default: Constant(-70.6 mV). Must return quantity with mV units when called with(shape,).I_ex_initializer (
Callable, optional) – Excitatory current initializer. Default: Constant(0.0 pA). Must return quantity with pA units when called with(shape,).I_in_initializer (
Callable, optional) – Inhibitory current initializer. Default: Constant(0.0 pA). Must return quantity with pA units when called with(shape,).w_initializer (
Callable, optional) – Adaptation current initializer. Default: Constant(0.0 pA). Must return quantity with pA units when called with(shape,).spk_fun (
Callable, optional) – Surrogate gradient function for spike generation. Default: ReluGrad(). Must be a differentiable spike function frombraintools.surrogate.spk_reset (
str, optional) – Spike reset mode. Default: ‘hard’. - ‘hard’: Stop gradient through reset (matches NEST behavior) - ‘soft’: Allow gradient through resetref_var (
bool, optional) – If True, expose boolean refractory state variable. Default: False. When True, createsself.refractoryindicating refractory status.name (
str, optional) – Model instance name. Default: None (auto-generated).
Parameter Mapping
Parameter
Default
Math equivalent
Description
in_size(required)
Population shape
V_peak0 mV
\(V_\mathrm{peak}\)
Spike detection threshold (if
Delta_T > 0)V_reset-60 mV
\(V_\mathrm{reset}\)
Reset potential
t_ref0 ms
\(t_\mathrm{ref}\)
Absolute refractory duration
g_L30 nS
\(g_\mathrm{L}\)
Leak conductance
C_m281 pF
\(C_\mathrm{m}\)
Membrane capacitance
E_L-70.6 mV
\(E_\mathrm{L}\)
Leak reversal potential
Delta_T2 mV
\(\Delta_T\)
Exponential slope factor
tau_w144 ms
\(\tau_w\)
Adaptation time constant
a4 nS
\(a\)
Subthreshold adaptation
b80.5 pA
\(b\)
Spike-triggered adaptation increment
V_th-50.4 mV
\(V_\mathrm{th}\)
Spike initiation threshold (in exponential term)
tau_syn_ex0.2 ms
\(\tau_{\mathrm{syn,ex}}\)
Excitatory exponential time constant
tau_syn_in2.0 ms
\(\tau_{\mathrm{syn,in}}\)
Inhibitory exponential time constant
I_e0 pA
\(I_\mathrm{e}\)
Constant external current
gsl_error_tol1e-6
(solver tolerance)
RKF45 local error tolerance
V_initializerConstant(-70.6 mV)
Membrane initializer
I_ex_initializerConstant(0 pA)
Excitatory current initializer
I_in_initializerConstant(0 pA)
Inhibitory current initializer
w_initializerConstant(0 pA)
Adaptation current initializer
spk_funReluGrad()
Surrogate spike function
spk_reset'hard'Reset mode; hard reset matches NEST behavior
ref_varFalseIf True, expose boolean refractory indicator
- V#
Membrane potential. Shape:
(*in_size,). Units: mV.- Type:
brainstate.HiddenState
- I_ex#
Excitatory synaptic current. Shape:
(*in_size,). Units: pA.- Type:
brainstate.HiddenState
- I_in#
Inhibitory synaptic current. Shape:
(*in_size,). Units: pA.- Type:
brainstate.HiddenState
- w#
Adaptation current. Shape:
(*in_size,). Units: pA.- Type:
brainstate.HiddenState
- refractory_step_count#
Remaining refractory steps. Shape:
(*in_size,). Dtype: int32.- Type:
brainstate.ShortTermState
- integration_step#
Persistent RKF45 internal step size. Shape:
(*in_size,). Units: ms.- Type:
brainstate.ShortTermState
- I_stim#
One-step delayed current buffer. Shape:
(*in_size,). Units: pA.- Type:
brainstate.ShortTermState
- last_spike_time#
Last emitted spike time. Shape:
(*in_size,). Units: ms. Updated tot + dton spike emission.- Type:
brainstate.ShortTermState
- refractory#
Boolean refractory indicator. Only exists if
ref_var=True. Shape:(*in_size,). Dtype: bool.- Type:
brainstate.ShortTermState, optional
- Raises:
If
V_reset >= V_peak- IfDelta_T < 0- IfV_peak < V_th- IfC_m <= 0- Ift_ref < 0- If any time constant<= 0- Ifgsl_error_tol <= 0- If(V_peak - V_th) / Delta_Tis too large (overflow risk in exponential term) - If numerical instability detected (V < -1e3or|w| > 1e6)
Notes
Implementation Details:
Adaptive integration: RKF45 adjusts step size \(h\) dynamically to meet error tolerance. Step size is persistent across simulation steps for efficiency.
Refractory semantics: During refractory, voltage is clamped to
V_resetin the ODE RHS anddV/dt = 0. This matches NEST exactly.Multiple spikes per step: With
t_ref = 0, multiple spikes can occur within one simulation step. Each spike triggers reset and adaptation increment.Overflow protection: Parameter validation checks that the exponential term \(\exp((V_{\text{peak}} - V_{th}) / \Delta_T)\) does not overflow.
Surrogate gradients: For backpropagation, spike generation uses
spk_fun(default: ReLU gradient). Hard reset (spk_reset='hard') stops gradient through reset, matching biological discontinuity.
Differences from other models:
aeif_cond_exp: Uses conductance-based synapses instead of current-based.aeif_psc_alpha: Uses alpha-function synapses instead of exponential.aeif_psc_delta: Uses delta-function (instantaneous) synapses.
Examples
Basic usage with constant input current:
>>> import brainpy.state as bp >>> import saiunit as u >>> import brainstate >>> >>> # Create population of 100 AdEx neurons >>> neurons = bp.aeif_psc_exp(100, I_e=200 * u.pA) >>> >>> # Initialize states >>> with brainstate.environ.context(dt=0.1 * u.ms): ... neurons.init_all_states() ... ... # Run for 100 ms ... spikes = [] ... for _ in range(1000): ... spike = neurons.update() ... spikes.append(spike)
With synaptic input and refractory period:
>>> # Create neurons with 2 ms refractory period >>> neurons = bp.aeif_psc_exp( ... in_size=100, ... t_ref=2.0 * u.ms, ... tau_syn_ex=5.0 * u.ms, ... tau_syn_in=10.0 * u.ms ... ) >>> >>> with brainstate.environ.context(dt=0.1 * u.ms): ... neurons.init_all_states() ... ... # Add excitatory input (positive weights) ... neurons.add_delta_input('exc', lambda: 100 * u.pA) ... ... # Simulation step ... spike = neurons.update(x=50 * u.pA) # External current
References
- get_spike(V=None)[source]#
Compute differentiable spike output using surrogate gradient.
Applies the surrogate gradient function
spk_funto a scaled voltage. The scaling maps the threshold region to a suitable range for the surrogate.- Parameters:
V (
ArrayLike, optional) – Membrane potential. Units: mV. If None, usesself.V.value. Shape:(*in_size,).- Returns:
Differentiable spike output in [0, 1]. Shape:
(*in_size,). Dtype: float. Values close to 1 indicate spike, close to 0 indicate silence.- Return type:
ArrayLike
Notes
The voltage is scaled by
(V - V_th) / (V_th - V_reset)before passing to the surrogate function. This normalizes the threshold crossing region for the surrogate gradient approximation.
- init_state(**kwargs)[source]#
Initialize persistent and short-term state variables.
- Parameters:
**kwargs – Unused compatibility parameters accepted by the base-state API.
- Raises:
ValueError – If an initializer cannot be broadcast to requested shape.
TypeError – If initializer outputs have incompatible units/dtypes for the corresponding state variables.
- update(x=Quantity(0., 'pA'))[source]#
Advance neuron state by one simulation timestep.
Performs one simulation step of the adaptive exponential integrate-and-fire neuron using adaptive RKF45 integration. Handles refractory clamping, spike detection, reset, adaptation increment, and synaptic input application.
The update sequence follows NEST semantics:
Integrate ODEs over \([t, t+dt]\) using adaptive RKF45:
Vectorized integration with adaptive step size
Inside integration: apply refractory clamp, detect spikes, reset voltage, increment adaptation, update refractory counter
Step size persists across simulation steps for efficiency
Post-integration processing:
Decrement refractory counter by 1
Apply delta inputs (spike weights) to \(I_{ex}\), \(I_{in}\)
Store external current \(x\) into one-step delayed buffer \(I_{\text{stim}}\)
Update
last_spike_timefor neurons that spiked
Return spike tensor:
Binary array indicating which neurons spiked during \([t, t+dt]\)
- Parameters:
x (
ArrayLike, optional) – External current input. Units: pA. Default: 0.0 pA. Shape: scalar or broadcastable to(*in_size,). Combined withcurrent_inputsand stored inI_stimfor next step.- Returns:
Binary spike tensor with dtype
jnp.float64and shapeself.V.value.shape. A value of1.0indicates at least one internal spike event occurred during the integrated interval \((t, t+dt]\).- Return type:
jax.Array- Raises:
ValueError – If numerical instability detected:
V < -1e3or|w| > 1e6.
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.See also
init_stateInitialize state variables
get_spikeCompute differentiable spike output