Blox Documentation

This page is a list of all of the Blox that Neuroblox can create.

A blox can be a neuron, a neural mass model, a receptor, a composite circuit of the preceding, a stimulus or current input, or experimental readouts (observers).

Neuron Blox

The following section documents some of the basic types of neurons that can be created, along with a description and parameters that can be set.

NeurobloxBasics.IFNeuronType
IFNeuron(name, namespace, C, θ, Eₘ, I_in)

Create a basic integrate-and-fire neuron. This follows Lapicque's equation (see Abbott [1], with parameters chosen to match the LIF/QIF neurons implemented as well):

\[\frac{dV}{dt} = \frac{I_{in} + jcn}{C}\]

where $jcn$ is any input to the blox.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • C: Membrane capicitance (μF).
  • θ: Threshold voltage (mV).
  • Eₘ: Resting membrane potential (mV).
  • I_in: External current input (μA).

References:

  1. Abbott, L. Lapicque's introduction of the integrate-and-fire model neuron (1907). Brain Res Bull 50, 303-304 (1999).
NeurobloxBasics.LIFNeuronType
LIFNeuron(name, namespace, C, θ, Eₘ, I_in)

Create a leaky integrate-and-fire neuron.

This largely follows the formalism and parameters given in Chapter 8 of Sterratt et al. [1], with the following equations:

\[\frac{dV}{dt} = \frac{\frac{-(V-E_m)}{R_m} + I_{in} + jcn}{C} \frac{dG}{dt} = -\frac{1}{\tau}G\]

where $jcn$ is any synaptic input to the blox (presumably a current G from another neuron).

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • C: Membrane capicitance (μF).
  • Eₘ: Resting membrane potential (mV).
  • Rₘ: Membrane resistance (kΩ).
  • τ: Synaptic time constant (ms).
  • θ: Threshold voltage (mV).
  • E_syn: Synaptic reversal potential (mV).
  • G_syn: Synaptic conductance (μA/mV).
  • I_in: External current input (μA).

References:

  1. Sterratt, D., Graham, B., Gillies, A., & Willshaw, D. (2011). Principles of Computational Modelling in Neuroscience. Cambridge University Press.
NeurobloxBasics.QIFNeuronType
QIFNeuron(; name, namespace, kwargs...)

Create a quadratic integrate-and-fire neuron.

The formal definition of this blox is:

\[ \frac{dV}{dt} = ((V - Eₘ)^2 / (Rₘ^2) + I_in + jcn) / C \frac{dG}{dt} = (-1 / τ₂)G + z \frac{dz}{dt} = (-1 / τ₁)z\]

where $jcn$ is any input to the blox.

Arguments:

  • C: Membrane capacitance (μF).
  • Rₘ: membrane resistance (kΩ).
  • E_syn: Synaptic reversal potential (mV).
  • G_syn: Synaptic conductance (mS).
  • τ₁: Timescale of decay of synaptic conductance (ms).
  • τ₂: Timescale of decay of synaptic spike variable (ms).
  • I_in: External current input (μA).
  • Eₘ: Resting membrane potential (mV).
  • Vᵣₑₛ: Post action potential (mV).
  • θ: Threshold potential (mV).
NeurobloxBasics.IzhikevichNeuronType
IzhikevichNeuron(; name, namespace, kwargs...)

Create an Izhikevich neuron, largely following the implementation in Chen and Campbell (2022), with synaptic decay.

The formal definition of this blox is:

\[\frac{dV}{dt} = V(V - α) - w + η + jcn \frac{dw}{dt} = a(bV - w) \frac{dG}{dt} = -(1 / τ)G + z \frac{dz}{dt} = -(1 / τ)z\]

Arguments:

  • α: The firing rate parameter (defaults to 0.6215).
  • η: Intrinsic current (defaults to 0.12 mA).
  • a: Timescale of the recovery variable (defaults to 0.0077 ms).
  • b: Sensitivity of the recovery variable to fluctuations in the voltage (defaults to −0.0062).
  • θ: Threshold voltage (defaults to 200.0 mV).
  • vᵣ: Reset potential after a spike (defaults to -200.0 mV).
  • wⱼ: The jump in the recovery variable after a spike (defaults to 0.0189).
  • sⱼ: Reset value for the synaptic spike variable after a spike (defaults to 1.2308).
  • gₛ: The synaptic conductance (defaults to 1.2308 mS).
  • eᵣ: The synaptic reversal potential (defaults to 1.0 mV).
  • τ: The synaptic decay time constant (defaults to 2.6 ms).

References:

  1. Izhikevich, E. (2003). Simple model of spiking neurons. IEEE Transactions on Neural Networks, 14(6), 1569–1572.
  2. Chen, L., & Campbell, S. A. (2022). Exact mean-field models for spiking neural networks with adaptation. Journal of Computational Neuroscience, 50(4), 445-469.
NeurobloxBasics.LIFExciNeuronType
LIFExciNeuron(; name, namespace, kwargs...)

Create an excitatory leaky integrate-and-fire neuron. This is a model that uses LIF equations for the voltage dynamics, but explicitly models the gating dynamics of glutamate and GABA receptors and their synaptic decay.

The formal definition of this blox is:

\[\frac{dV}{dt} = (-g_L(V - V_L) - S_\text{AMPA, ext} g_\text{AMPA, ext} (V - V_E) - S_\text{GABA} g_\text{GABA}(V - V_I) - S_\text{AMPA} g_\text{AMPA} (V - V_E) - jcn) / C D(S_\text{AMPA}) = - S_\text{AMPA} / τ_\text{AMPA}\\ D(S_\text{GABA}) = - S_\text{GABA} / τ_\text{GABA}\\ D(S_\text{NMDA}) = - S_\text{NMDA} / τ_\text{NMDA, decay} + αx(1 - S_\text{NMDA}) \frac{dx}{dt} = -x / τ_\text{NMDA, rise} \frac{dS_\text{AMPA, ext}}{dt} = -S_\text{AMPA, ext} / τ_\text{AMPA}\]

Arguments:

  • g_L: Leak conductance (mS).
  • V_L: Leak reversal potential (mV).
  • V_E: Excitatory reversal potential (mV).
  • V_I: Inhibitory reversal potential (mV).
  • θ: Firing threshold potential (mV).
  • V_reset: Reset potential after a spike (mV).
  • C: Membrane capacitance (μF).
  • τ_AMPA: time constant for the closing of AMPA receptors (ms).
  • τ_GABA: time constant for the closing of GABA receptors (ms).
  • τNMDAdecay: time constant for the closing of NMDA receptors (ms).
  • τNMDArise: time constant for the decay of NMDA current post spike (ms).
  • t_refract: Refractory period after a spike (ms).
  • α: Scaling for the rise of NMDA current (ms^-1).
  • g_AMPA: Synaptic conductance for AMPA glutamate receptors (mS).
  • gAMPAext: Synaptic conductance for external current input through AMPA receptors (mS).
  • g_GABA: Synaptic conductance for GABA glutamate receptors (mS).
  • g_NMDA: Synaptic conductance for NMDA glutamate receptors (mS).
  • Mg: Magnesium ion concentration (mM).
  • exciscalingfactor: Excitatory scaling factor.
  • inhscalingfactor: Inhibitory scaling factor.
NeurobloxBasics.LIFInhNeuronType
LIFInhNeuron(; name, namespace, kwargs...)

Create an inhibitory leaky integrate-and-fire neuron. This is a model that uses LIF equations for the voltage dynamics, but explicitly models the gating dynamics of AMPA and GABA receptors and their synaptic decay.

The formal definition of this blox is:

\[\frac{dV}{dt} = (1 - \mathbb{1}_\text{refrac}) (-g_L(V - V_L) - S_\text{AMPA, ext} g_\text{AMPA, ext} (V - V_E) - S_\text{GABA} g_\text{GABA}(V - V_I) - S_\text{AMPA} g_\text{AMPA} (V - V_E) - jcn) / C D(S_\text{AMPA}) = - S_\text{AMPA} / τ_\text{AMPA}\ D(S_\text{GABA}) = - S_\text{GABA} / τ_\text{GABA}\ D(S_\text{AMPA, ext}) = - S_\text{AMPA, ext} / τ_\text{AMPA}\]

Arguments:

  • g_L: Leak conductance (mS).
  • V_L: Leak reversal potential (mV).
  • V_E: Excitatory reversal potential (mV).
  • V_I: Inhibitory reversal potential (mV).
  • θ: Firing threshold potential (mV).
  • V_reset: Reset potential after a spike (mV).
  • C: Membrane capacitance (μF).
  • τ_AMPA: time constant for the closing of AMPA receptors (ms).
  • τ_GABA: time constant for the closing of GABA receptors (ms).
  • t_refract: Refractory period after a spike (ms).
  • α: Scaling for the rise of NMDA current (ms^-1).
  • g_AMPA: Synaptic conductance for AMPA glutamate receptors (mS).
  • gAMPAext: Synaptic conductance for external current input through AMPA receptors (mS).
  • g_GABA: Synaptic conductance for GABA glutamate receptors (mS).
  • g_NMDA: Synaptic conductance for NMDA glutamate receptors (mS).
  • Mg: Magnesium ion concentration (mM).
  • exciscalingfactor: Excitatory scaling factor.
  • inhscalingfactor: Inhibitory scaling factor.
NeurobloxBasics.PINGNeuronExciType
PINGNeuronExci(name, namespace, C, g_Na, V_Na, g_K, V_K, g_L, V_L, I_ext, τ_R, τ_D)

Create an excitatory neuron from Borgers et al. (2008). The formal definition of this blox is:

\[\frac{dV}{dt} = \frac{1}{C}(-g_{Na}*m_{\infty}^3*h*(V - V_{Na}) - g_K*n^4*(V - V_K) - g_L*(V - V_L) + I_{ext} + jcn) \m_{\infty} = \frac{a_m(V)}{a_m(V) + b_m(V)} \frac{dn}{dt} = a_n(V)*(1 - n) - b_n(V)*n \frac{dh}{dt} = a_h(V)*(1 - h) - b_h(V)*h \frac{ds}{dt} = \frac{1}{2}*(1 + \tanh(V/10))*(\frac{1 - s}{\tau_R} - \frac{s}{\tau_D})\]

where $jcn$ is any input to the blox. Note that this is a modified Hodgkin-Huxley formalism with an additional synaptic accumulation term. Synapses are added into the $jcn$ term by connecting the postsynaptic neuron's voltage to the presynaptic neuron's output:

\[jcn = w*s*(V_E - V)\]

where $w$ is the weight of the synapse and $V_E$ is the reversal potential of the excitatory synapse.

Inputs:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • C: Membrane capacitance (defaults to 1.0).
  • g_Na: Sodium conductance (defaults to 100.0).
  • V_Na: Sodium reversal potential (defaults to 50.0).
  • g_K: Potassium conductance (defaults to 80.0).
  • V_K: Potassium reversal potential (defaults to -100.0).
  • g_L: Leak conductance (defaults to 0.1).
  • V_L: Leak reversal potential (defaults to -67.0).
  • I_ext: External current (defaults to 0.0).
  • τ_R: Rise time of synaptic conductance (defaults to 0.2).
  • τ_D: Decay time of synaptic conductance (defaults to 2.0).
NeurobloxBasics.PINGNeuronInhibType
PINGNeuronInhib(name, namespace, C, g_Na, V_Na, g_K, V_K, g_L, V_L, I_ext, τ_R, τ_D)

Create an inhibitory neuron from Borgers et al. (2008). The formal definition of this blox is:

\[\frac{dV}{dt} = \frac{1}{C}(-g_{Na}*m_{\infty}^3*h*(V - V_{Na}) - g_K*n^4*(V - V_K) - g_L*(V - V_L) + I_{ext} + jcn) \m_{\infty} = \frac{a_m(V)}{a_m(V) + b_m(V)} \frac{dn}{dt} = a_n(V)*(1 - n) - b_n(V)*n \frac{dh}{dt} = a_h(V)*(1 - h) - b_h(V)*h \frac{ds}{dt} = \frac{1}{2}*(1 + \tanh(V/10))*(\frac{1 - s}{\tau_R} - \frac{s}{\tau_D})\]

where $jcn$ is any input to the blox. Note that this is a modified Hodgkin-Huxley formalism with an additional synaptic accumulation term. Synapses are added into the $jcn$ term by connecting the postsynaptic neuron's voltage to the presynaptic neuron's output:

\[jcn = w*s*(V_I - V)\]

where $w$ is the weight of the synapse and $V_I$ is the reversal potential of the inhibitory synapse.

Inputs:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • C: Membrane capacitance (defaults to 1.0).
  • g_Na: Sodium conductance (defaults to 35.0).
  • V_Na: Sodium reversal potential (defaults to 55.0).
  • g_K: Potassium conductance (defaults to 9.0).
  • V_K: Potassium reversal potential (defaults to -90.0).
  • g_L: Leak conductance (defaults to 0.1).
  • V_L: Leak reversal potential (defaults to -65.0).
  • I_ext: External current (defaults to 0.0).
  • τ_R: Rise time of synaptic conductance (defaults to 0.5).
  • τ_D: Decay time of synaptic conductance (defaults to 10.0).
Missing docstring.

Missing docstring for HHNeuronExci. Check Documenter's build log for details.

Missing docstring.

Missing docstring for HHNeuronInhib. Check Documenter's build log for details.

Missing docstring.

Missing docstring for HHNeuronFSI. Check Documenter's build log for details.

NeurobloxPharma.MetabolicHHNeuronType
MetabolicHHNeuron(name, namespace, neurontype,
	Naᵢᵧ, ρₘₐₓ, α, λ, ϵ₀, O₂ᵦ, γ, β, ϵₖ, Kₒᵦ, Gᵧ, Clᵢ, Clₒ, R, T, F,
	Gₙₐ, Gₖ, Gₙₐ_L, Gₖ_L, G_cl_L, C_m, I_in, G_exc, G_inh, E_syn_exc, E_syn_inh)

Create a Metabolic Hodgkin-Huxley Neuron. This model accounts for
dynamic ion concentrations, oxygen consumption and astrocytic buffering.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • neurontype: excitatory or inhibitory.
  • Naᵢᵧ: Intracellular Na+ concentration (mM).
  • ρₘₐₓ: Maximum pump rate (mM/s).
  • α: Conversion factor from pump current to O2 consumption rate (g/mol).
  • λ: Relative cell density.
  • ϵ₀: O2 diffusion rate (s^-1).
  • O₂ᵦ: O2 buffer concentration (mg/L).
  • γ: Conversion factor from current to concentration (mM/s)/(uA/cm2).
  • β: Ratio of intracellular vs extracellular volume.
  • ϵₖ: K+ diffusion rate (1/s).
  • Kₒᵦ: K+ buffer concentration (mM).
  • Gᵧ: Glia uptake strength of K+ (mM/s).
  • Clᵢ: Intracellular Cl- concentration (mM).
  • Clₒ: Extracellular Cl- concentration (mM).
  • R: Ideal gas constant (J/(mol*K)).
  • T: Temperature (K).
  • F: Faraday's constant (C/mol).
  • Gₙₐ: Na+ maximum conductance (mS/cm^2).
  • Gₖ: K+ maximum conductance (mS/cm^2).
  • Gₙₐ_L: Na+ leak conductance (mS/cm^2).
  • Gₖ_L: K+ leak conductance (mS/cm^2).
  • GclL: Cl- leak conductance (mS/cm^2).
  • C_m: Membrane capacitance (uF/cm^2).
  • I_in: External current input (uA/cm^2).
  • G_exc: Conductance of excitatory synapses (mS/cm^2).
  • G_inh: Conductance of inhibitory synapses (mS/cm^2).
  • Esynexc: Excitatory synaptic reversal potential (mV).
  • Esyninh: Inhibitory synaptic reversal potential (mV).

References:

  1. Dutta, Shrey, et al. "Mechanisms underlying pathological cortical bursts during metabolic depletion." Nature Communications 14.1 (2023): 4792.

Neural Mass Model Blox

Neural mass models simulate the average activity and voltage of a collection of neurons.

NeurobloxBasics.NGNMM_IzhType
NGNMM_Izh(name, namespace, ...)

This is the basic Izhikevich next-gen neural mass as described in [1]. The corresponding connector is set up to allow for connections between masses, but the user must add their own $ \kappa $ values to the connection weight as there is no good way of accounting for this weight within/between regions.

Currently, the connection weights include the presynaptic $ g_s $, but this could be changed.

Equations: To be added once we have a final form that we like here.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • Other parameters: See reference for full list. Note that parameters are scaled so that units of time are in milliseconds.

Citation:

  1. Chen, L., & Campbell, S. A. (2022). Exact mean-field models for spiking neural networks with adaptation. Journal of Computational Neuroscience, 50(4), 445-469.
NeurobloxBasics.NGNMM_QIFType
NGNMM_QIF(name, namespace, ...)

This is the basic QIF next-gen neural mass as described in [1]. This includes the connections via firing rate as described in [1] and the optional noise term.

Equations: To be added once we have a final form that we like here.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • Other parameters: See reference for full list. Note that parameters are scaled so that units of time are in milliseconds.

Citation: Theta-nested gamma bursts by Torcini group.

NeurobloxBasics.LinearNeuralMassType
LinearNeuralMass(name, namespace)

Create standard linear neural mass blox with a single internal state. There are no parameters in this blox. This is a blox of the sort used for spectral DCM modeling. The formal definition of this blox is:

\[\frac{d}{dx} = \sum{jcn}\]

where $jcn$` is any input to the blox.

Arguments:

  • name: Options containing specification about deterministic.
  • namespace: Additional namespace above name if needed for inheritance.
NeurobloxBasics.HarmonicOscillatorType
HarmonicOscillator(name, namespace, ω, ζ, k, h)

Create a harmonic oscillator blox with the specified parameters.

The formal definition of this blox is:

\[\frac{dx}{dt} = y-(2*\omega*\zeta*x)+ k*(2/\pi)*(atan((\sum{jcn})/h) \frac{dy}{dt} = -(\omega^2)*x\]

where $jcn$ is any input to the blox.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • ω: Base frequency. Note the default value is scaled to give oscillations in milliseconds to match other blocks.
  • ζ: Damping ratio.
  • k: Gain.
  • h: Threshold.
NeurobloxBasics.JansenRitType
JansenRit(name, namespace, τ, H, λ, r, cortical, delayed)

Create a Jansen Rit blox as described in Liu et al. The formal definition of this blox is:

\[\frac{dx}{dt} = y-\frac{2}{\tau}x \frac{dy}{dt} = -\frac{x}{\tau^2} + \frac{H}{\tau} [\frac{2\lambda}{1+\text{exp}(-r*\sum{jcn})} - \lambda]\]

where $jcn$ is any input to the blox.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • τ: Time constant. Defaults to 1 for cortical regions, 14 for subcortical.
  • H: See equation for use. Defaults to 0.02 for both cortical and subcortical regions.
  • λ: See equation for use. Defaults to 5 for cortical regions, 400 for subcortical.
  • r: See equation for use. Defaults to 0.15 for cortical regions, 0.1 for subcortical.
  • cortical: Boolean to determine whether to use cortical or subcortical parameters. Specifying any of the parameters above will override this.
  • delayed: Boolean to indicate whether states are delayed

Citations:

  1. Liu C, Zhou C, Wang J, Fietkiewicz C, Loparo KA. The role of coupling connections in a model of the cortico-basal ganglia-thalamocortical neural loop for the generation of beta oscillations. Neural Netw. 2020 Mar;123:381-392. doi: 10.1016/j.neunet.2019.12.021.
NeurobloxBasics.WilsonCowanType
WilsonCowan(name, namespace, τ_E, τ_I, a_E, a_I, c_EE, c_IE, c_EI, c_II, θ_E, θ_I, η)

Create a standard Wilson Cowan blox. The formal definition of this blox is:

\[\frac{dE}{dt} = \frac{-E}{\tau_E} + \frac{1}{1 + \text{exp}(-a_E*(c_{EE}*E - c_{IE}*I - \theta_E + \eta*(\sum{jcn}))} \frac{dI}{dt} = \frac{-I}{\tau_I} + \frac{1}{1 + exp(-a_I*(c_{EI}*E - c_{II}*I - \theta_I)}\]

where $jcn$ is any input to the blox.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • Others: See equation for use.
NeurobloxBasics.LarterBreakspearType
LarterBreakspear(name, namespace, ...)

Create a Larter Breakspear blox described in Endo et al. For a full list of the parameters used see the reference. If you need to modify the parameters, see Chesebro et al. and van Nieuwenhuizen et al. for physiological ranges.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • Other parameters: See reference for full list. Note that parameters are scaled so that units of time are in milliseconds.

Citations:

  1. Endo H, Hiroe N, Yamashita O. Evaluation of Resting Spatio-Temporal Dynamics of a Neural Mass Model Using Resting fMRI Connectivity and EEG Microstates. Front Comput Neurosci. 2020 Jan 17;13:91. doi: 10.3389/fncom.2019.00091.
  2. Chesebro AG, Mujica-Parodi LR, Weistuch C. Ion gradient-driven bifurcations of a multi-scale neuronal model. Chaos Solitons Fractals. 2023 Feb;167:113120. doi: 10.1016/j.chaos.2023.113120.
  3. van Nieuwenhuizen, H, Chesebro, AG, Polis, C, Clarke, K, Strey, HH, Weistuch, C, Mujica-Parodi, LR. Ketosis regulates K+ ion channels, strengthening brain-wide signaling disrupted by age. Preprint. bioRxiv 2023.05.10.540257; doi: https://doi.org/10.1101/2023.05.10.540257.
NeurobloxBasics.KuramotoOscillatorType
KuramotoOscillator(name, namespace, ...)

Simple implementation of the Kuramoto oscillator as described in the original paper [1]. Useful for general models of synchronization and oscillatory behavior. The general form of the Kuramoto oscillator is given by: Equations:

\[ \begin{equation} \dot{\theta_i} = \omega_i + \frac{1}{N}\sum_{j=1}^N{K_{i, j}\text{sin}(\theta_j - \theta_i)} \end{equation}\]

Where this describes the connection between regions $i$ and $j$. An alternative form which includes a noise term for each region is also provided, taking the form:

\[ \begin{equation} \dot{\theta_i} = \omega_i + \zeta dW_i \frac{1}{N}\sum_{j=1}^N{K_{i, j}\text{sin}(\theta_j - \theta_i)} \end{equation}\]

where $W_i$ is a Wiener process and $\zeta_i$ is the noise strength.

Keyword arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • include_noise (default false) determines if brownian noise is included in the dynamics of the blox.
  • Other parameters: See reference for full list. Note that parameters are scaled so that units of time are in milliseconds. Default parameter values are taken from [2].

Citations:

  1. Kuramoto, Y. (1975). Self-entrainment of a population of coupled non-linear oscillators. In: Araki, H. (eds) International Symposium on Mathematical Problems in Theoretical Physics. Lecture Notes in Physics, vol 39. Springer, Berlin, Heidelberg. https://doi.org/10.1007/BFb0013365

  2. Sermon JJ, Wiest C, Tan H, Denison T, Duchet B. Evoked resonant neural activity long-term dynamics can be reproduced by a computational model with vesicle depletion. Neurobiol Dis. 2024 Jun 14;199:106565. doi: 10.1016/j.nbd.2024.106565. Epub ahead of print. PMID: 38880431.

NeurobloxBasics.VanDerPolType
VanDerPol(; name, namespace = nothing, θ=1.0, ϕ=0.1, include_noise = false)

Create a neural mass model whose activity variable follows the dynamics of the (stochastic) van der Pol oscillator.

The formal definition of this blox is:

\[\frac{dx}{dt} = y \frac{dy}{dt} = θ(1 - x^2)y - x + ϕξ + jcn\]

where jcn is any input to the blox, and ξ is a Brownian variable that is added if include_noise is set.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • θ: damping strength
  • ϕ: strength of the Brownian motion
NeurobloxBasics.Generic2dOscillatorType
Generic2dOscillator(name, namespace, ...)

The Generic2dOscillator model is a generic dynamic system with two state variables. The dynamic equations of this model are composed of two ordinary differential equations comprising two nullclines. The first nullcline is a cubic function as it is found in most neuron and population models; the second nullcline is arbitrarily configurable as a polynomial function up to second order. The manipulation of the latter nullcline's parameters allows to generate a wide range of different behaviours.

Equations:

\[ \begin{align} \dot{V} &= d \, \tau (-f V^3 + e V^2 + g V + \alpha W + \gamma I) \\ \dot{W} &= \dfrac{d}{ au}\,\,(c V^2 + b V - \beta W + a) \end{align}\]

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • Other parameters: See reference for full list. Note that parameters are scaled so that units of time are in milliseconds.

Citations: FitzHugh, R., Impulses and physiological states in theoretical models of nerve membrane, Biophysical Journal 1: 445, 1961.

Nagumo et.al, An Active Pulse Transmission Line Simulating Nerve Axon, Proceedings of the IRE 50: 2061, 1962.

Stefanescu, R., Jirsa, V.K. Reduced representations of heterogeneous mixed neural networks with synaptic coupling. Physical Review E, 83, 2011.

Jirsa VK, Stefanescu R. Neural population modes capture biologically realistic large-scale network dynamics. Bulletin of Mathematical Biology, 2010.

Stefanescu, R., Jirsa, V.K. A low dimensional description of globally coupled heterogeneous neural networks of excitatory and inhibitory neurons. PLoS Computational Biology, 4(11), 2008).

NeurobloxBasics.OUProcessType
OUProcess(; name, namespace, μ, σ, τ)

Create a neural mass model whose activity variable follows an Ornstein-Uhlenbeck process.

The formal definition of this blox is:

\[\frac{dx}{dt} = (-x + μ + jcn)/τ + \sqrt{2 / τ} σw\]

Where w is a Brownian variable and jcn is the input to the blox.

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • τ: relaxation time
  • μ: Mean of the OU process
  • σ: Strength of the Brownian motion (variance of OUProcess process is τ*σ^2/2)
NeurobloxPharma.NGNMM_thetaType
NGNMM_theta(name, namespace, ...)
    Create a next-gen neural mass model of coupled theta neuron populations. For a full list of the parameters used see the reference.
    Each mass consists of a population of two neurons ``a`` and ``b``, coupled using different synaptic terms ``g``. The entire expression of these is given by:

\[ \frac{a_e}{dt} = \frac{1}{C_e}(b_e*(a_e-1) - (\Delta_e/2)*((a_e+1)^2-b_e^2) - \eta_{0e}*b_e*(a_e+1) - (v_{syn, ee}*g_{ee}+v_{syn, ei}*g_{ei})*(b_e*(a_e+1)) - (g_{ee}/2+g_{ei}/2)*(a_e^2-b_e^2-1)) \frac{b_e}{dt} = \frac{1}{C_e}*((b_e^2-(a_e-1)^2)/2 - \Delta_e*b_e*(a_e+1) + (\eta_{0e}/2)*((a_e+1)^2-b_e^2) + (v_{syn, ee}*(g_{ee}/2)+v_{syn, ei}*(g_{ei}/2))*((a_e+1)^2-b_e^2) - a_e*b_e*(g_{ee}+g_{ei})) \frac{a_i}{dt} = \frac{1}{C_i}(b_i*(a_i-1) - (\Delta_i/2)*((a_i+1)^2-b_i^2) - \eta_{0i}*b_i*(a_i+1) - (v_{syn, ie}*g_{ie}+v_{syn, ii}*g_{ii})*(b_i*(a_i+1)) - (g_{ie}/2+g_{ii}/2)*(a_i^2-b_i^2-1)) \frac{b_i}{dt} = \frac{1}{C_i}*((b_i^2-(a_i-1)^2)/2 - \Delta_i*b_i*(a_i+1) + (\eta_{0i}/2)*((a_i+1)^2-b_i^2) + (v_{syn, ie}*(g_{ie}/2)+v_{syn, ii}*(g_{ii}/2))*((a_i+1)^2-b_i^2) - a_i*b_i*(g_{ie}+g_{ii})) \frac{g_ee}{dt} = \alpha_{inv, ee} (\frac{k_{ee}}{C_e \pi} \frac{1-a_e^2-b_e^2}{(1+2*a_e+a_e^2+b_e^2)} - g_{ee}) \frac{g_ei}{dt} = \alpha_{inv, ei} (\frac{k_{ei}}{C_i \pi} \frac{1-a_i^2-b_i^2}{(1+2*a_i+a_i^2+b_i^2)} - g_{ei}) \frac{g_ie}{dt} = \alpha_{inv, ie} (\frac{k_{ie}}{C_e \pi} \frac{1-a_e^2-b_e^2}{(1+2*a_e+a_e^2+b_e^2)} - g_{ie}) \frac{g_ii}{dt} = \alpha_{inv, ii} (\frac{k_{ii}}{C_i \pi} \frac{1-a_i^2-b_i^2}{(1+2*a_i+a_i^2+b_i^2)} - g_{ii})\]

Can alternatively be called by $NextGenerationEI()$, but this is deprecated and will be removed in future updates.

Citations:

  1. Byrne Á, O'Dea RD, Forrester M, Ross J, Coombes S. Next-generation neural mass and field modeling. J Neurophysiol. 2020 Feb 1;123(2):726-742. doi: 10.1152/jn.00406.2019.
Missing docstring.

Missing docstring for NextGenerationEI. Check Documenter's build log for details.

Composite Blox

Composite blox typically represent larger-scale brain structures.

Pharmaceutical Module

The following composite blox are used for simulations of the effects of drug dosing on neural circuits.

NeurobloxPharma.WinnerTakeAllType
WinnerTakeAll

Creates a winner-take-all local circuit found in neocortex, typically 5 pyramidal (excitatory) neurons send synapses to a single interneuron (inhibitory) and receive feedback inhibition from that interneuron.

References:

  • Coultrip, Robert, Richard Granger, and Gary Lynch. “A Cortical Model of Winner-Take-All Competition via Lateral Inhibition.” Neural Networks 5, no. 1 (January 1, 1992): 47-54.
  • Pathak A., Brincat S., Organtzidis H., Strey H., Senneff S., Antzoulatos E., Mujica-Parodi L., Miller E., Granger R. "Biomimetic model of corticostriatal micro-assemblies discovers new neural code.", bioRxiv 2023.11.06.565902, 2024
Missing docstring.

Missing docstring for Cortical. Check Documenter's build log for details.

Missing docstring.

Missing docstring for LateralAmygdalaCluster. Check Documenter's build log for details.

Missing docstring.

Missing docstring for LateralAmygdala. Check Documenter's build log for details.

Missing docstring.

Missing docstring for CentralAmygdala. Check Documenter's build log for details.

Missing docstring.

Missing docstring for Striatum. Check Documenter's build log for details.

Missing docstring.

Missing docstring for GPi. Check Documenter's build log for details.

Missing docstring.

Missing docstring for GPe. Check Documenter's build log for details.

Missing docstring.

Missing docstring for Thalamus. Check Documenter's build log for details.

Missing docstring.

Missing docstring for STN. Check Documenter's build log for details.

Missing docstring.

Missing docstring for Nuc_Reticularis. Check Documenter's build log for details.

Deep Brain Stimulation

The following composite blox are used for simulations of deep brain stimulation.

Missing docstring.

Missing docstring for HHNeuronInhib_MSN_Adam. Check Documenter's build log for details.

Missing docstring.

Missing docstring for HHNeuronInhib_FSI_Adam. Check Documenter's build log for details.

Missing docstring.

Missing docstring for HHNeuronExci_STN_Adam. Check Documenter's build log for details.

Missing docstring.

Missing docstring for HHNeuronInhib_GPe_Adam. Check Documenter's build log for details.

Missing docstring.

Missing docstring for Striatum_MSN_Adam. Check Documenter's build log for details.

Missing docstring.

Missing docstring for Striatum_FSI_Adam. Check Documenter's build log for details.

Missing docstring.

Missing docstring for GPe_Adam. Check Documenter's build log for details.

Missing docstring.

Missing docstring for STN_Adam. Check Documenter's build log for details.

Discrete Blox

Discrete blox are ones in which the state changes discretely, rather than according to a differential equation.

Missing docstring.

Missing docstring for Matrisome. Check Documenter's build log for details.

Missing docstring.

Missing docstring for Striosome. Check Documenter's build log for details.

Missing docstring.

Missing docstring for TAN. Check Documenter's build log for details.

Missing docstring.

Missing docstring for SNc. Check Documenter's build log for details.

Receptor Blox

Missing docstring.

Missing docstring for MoradiNMDAR. Check Documenter's build log for details.

Source Blox

Sources are external inputs that feed into neural circuits. These may represent electrical inputs or sensory stimuli.

External Current Inputs

NeurobloxBasics.ConstantInputType
ConstantInput(; namespace, namespace, I)

Create a constant input current.

Inputs:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • I: External current input (mA).
NeurobloxBasics.PoissonSpikeTrainType
PoissonSpikeTrain(rate, tspan; name, namespace, N_trains, prob_dt, rng)

Create an input that generates spikes according to a Poisson process.

Inputs:

  • rate: A single rate or vector of rates for the Poisson processes.
  • tspan: A single timespan tuple or vector of timespan tuples.
  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • N_trains: Number of generated spike trains.
  • prob_dt: Timestep of the problem, used to determine the time step for the Poisson process.
  • rng: Random number generator.
Missing docstring.

Missing docstring for VoltageClampSource. Check Documenter's build log for details.

Missing docstring.

Missing docstring for PulsesInput. Check Documenter's build log for details.

Deep Brain Stimulation Stimuli

NeurobloxDBS.DBSType
DBS(; name, namespace=nothing, frequency=130.0, amplitude=2.5, pulse_width=0.066, 
    offset=0.0, start_time=0.0, smooth=1e-4)

Create a continuous deep brain stimulation (DBS) stimulus with regular pulses.

Arguments:

  • name: Name given to ODESystem object within the blox
  • namespace: Additional namespace above name if needed for inheritance
  • frequency: Pulse frequency in Hz
  • amplitude: Pulse amplitude in arbitrary units
  • pulse_width: Duration of each pulse in ms
  • offset: Baseline value of the signal between pulses
  • start_time: Time delay before stimulation begins in ms
  • smooth: Smoothing parameter for pulse transitions, set to 0 for sharp transitions

Returns a DBS stimulus blox that outputs square pulses with specified parameters.

NeurobloxDBS.ProtocolDBSType
ProtocolDBS(; name, namespace=nothing, frequency=130.0, amplitude=2.5,
              pulse_width=0.066, offset=0.0, start_time=0.0, smooth=1e-4,
              pulses_per_burst=10, bursts_per_block=12, 
              pre_block_time=200.0, inter_burst_time=200.0)

Create a deep brain stimulation (DBS) stimulus consisting of a block of pulse bursts.

Arguments:

  • name: Name given to ODESystem object within the blox
  • namespace: Additional namespace above name if needed for inheritance
  • frequency: Pulse frequency in Hz
  • amplitude: Pulse amplitude in arbitrary units
  • pulse_width: Duration of each pulse in ms
  • offset: Baseline value of the signal between pulses
  • start_time: Time delay before stimulation begins in ms
  • smooth: Smoothing parameter for pulse transitions, set to 0 for sharp transitions
  • pulsesperburst: Number of pulses in each burst
  • burstsperblock: Number of bursts in the block
  • preblocktime: Time before the block starts in ms
  • interbursttime: Time between bursts in ms

Returns a DBS stimulus blox that outputs a block of pulse bursts.

Missing docstring.

Missing docstring for SquareStimulus. Check Documenter's build log for details.

Missing docstring.

Missing docstring for BurstStimulus. Check Documenter's build log for details.

Sensory Stimuli

Missing docstring.

Missing docstring for ImageStimulus. Check Documenter's build log for details.

Observers

Observers are experimental readouts from neural simulations, such as fMRI signals.

NeurobloxBasics.BalloonModelType
BalloonModel(;name, namespace=nothing, lnκ=0.0, lnτ=0.0)

Create a balloon model blox which computes the hemodynamic responses to some underlying neuronal activity
The formal definition of this blox is:
```math
    \frac{ds}{dt} = \text{jcn} - \kappa s - \gamma (u-1)
    \frac{du}{dt} = s
    \frac{d\nu}{dt} = u - v^{1/\alpha}
    \frac{dq}{dt} = u E(u, E_0)/E_0 - v^{1/\alpha} q/v
```
where ``jcn`` is any input to the blox (represents the neuronal activity).

Includes also the bold measurement function:
```math
    \lambda(\nu, q) = V_0 \left[ k_1 (1-q) + k_2 \left( 1 - \frac{q}{v} \right) + k_3 (1-v)\right]
```

Arguments:

  • name: Name given to ODESystem object within the blox.
  • namespace: Additional namespace above name if needed for inheritance.
  • lnκ: logarithmic prefactor to signal decay H[1], set to 0 for standard parameter value.
  • lnτ: logarithmic prefactor to transit time H[3], set to 0 for standard parameter value.
  • lnϵ: logarithm of ratio of intra- to extra-vascular signal

NB: the prefix ln of the variables u, ν, q as well as the parameters κ, τ denotes their transformation into logarithmic space to enforce their positivity. This transformation is considered in the derivates of the model equations below.

Citations:

  1. Stephan K E, Weiskopf N, Drysdale P M, Robinson P A, and Friston K J. Comparing Hemodynamic Models with DCM. NeuroImage 38, no. 3 (2007): 387–401. doi: 10.1016/j.neuroimage.2007.07.040
  2. Hofmann D, Chesebro A G, Rackauckas C, Mujica-Parodi L R, Friston K J, Edelman A, and Strey H H. Leveraging Julia's Automated Differentiation and Symbolic Computation to Increase Spectral DCM Flexibility and Speed, 2023. doi: 10.1101/2023.10.27.564407
NeurobloxBasics.boldsignal_endo_balloonFunction

Simulated BOLD response to input

This is the version of the Balloon-Windkessel model used in Endo et al. (2020), and subsequently in Antal et al. (2023) and van Niuewenhuizen et al. (2023).
As it is only a convolution after the simulation of the neural activity, it is not a block for inclusion in the ODESystem.