Extracellular reaction-diffusion tutorial
We have expanded the capabilities the NEURON reaction diffusion module to support a macroscopic model of the extracellular space, described in a recent paper. Here is brief a tutorial that provides an overview of the python interface.
Download files
This tutorial uses two files. Download them and save them to your working directory:
The first is simple mechanism written for this tutorial that releases potassium at a constant rate per surface area. The second is a morphology from NeuroMorpho.Org.
Begin by compiling the mod file. On Linux and Mac, this can be done by running nrnivmodl from the command line. The same works for Windows beginning with NEURON 7.7; alternatively, on Windows regardless of NEURON version, one can use the graphical tool mknrndll.
Getting started
from neuron import h, rxd
(If you're using NEURON 7.6.x, use "from neuron import h, crxd as rxd" instead. The above is the recommended usage for 7.7+.)
For reaction diffusion in the extracellular space (ECS) you have to import the crxd from neuron. The first example we will place two single compartment neurons in a closed box of extracellular space.
rxd.options.enable.extracellular=True
To use extracellular reaction diffusion it must be enabled, it is enabled by default in NEURON 7.7 or later.
rxd.nthread(4)
Extracellular rxd supports multithreaded parallelization.
h.load_file('stdrun.hoc')
h.load_file('import3d.hoc')
class Cell:
def __init__(self,filename):
"""Read geometry from a given SWC file and create a cell with a K+ source"""
cell = h.Import3d_SWC_read()
cell.input(filename)
h.Import3d_GUI(cell, 0)
i3d = h.Import3d_GUI(cell, 0)
i3d.instantiate(self)
for sec in self.all:
sec.nseg = 1 + 10 * int(sec.L / 5)
sec.insert('steady_k')
def extrema(self):
"""Give the bounding box that contains the cell"""
xlo = ylo = zlo = xhi = yhi = zhi = None
for sec in self.all:
n3d = sec.n3d()
xs = [sec.x3d(i) for i in range(n3d)]
ys = [sec.y3d(i) for i in range(n3d)]
zs = [sec.z3d(i) for i in range(n3d)]
my_xlo, my_ylo, my_zlo = min(xs), min(ys), min(zs)
my_xhi, my_yhi, my_zhi = max(xs), max(ys), max(zs)
if xlo is None:
xlo, ylo, zlo = my_xlo, my_ylo, my_zlo
xhi, yhi, zhi = my_xhi, my_yhi, my_zhi
else:
xlo, ylo, zlo = min(xlo, my_xlo), min(ylo, my_ylo), min(zlo, my_zlo)
xhi, yhi, zhi = max(xhi, my_xhi), max(yhi, my_yhi), max(zhi, my_zhi)
return (xlo, ylo, zlo, xhi, yhi, zhi)
mycell = Cell('c91662.CNG.swc')
Create a cell from the given SWC file, in this example we use a CA1 pyramidal cell taken from NeuroMorpho.Org. We insert a current source steady_k
into all sections as a simple proof of concept example.
xlo, ylo, zlo, xhi, yhi, zhi = mycell.extrema()
padding = 50
ecs = rxd.Extracellular(xlo - padding, ylo - padding, zlo - padding,
xhi + padding, yhi + padding, zhi + padding,
dx=(20, 20, 50), volume_fraction=0.2, tortuosity=1.6)
We define the extracellular region by specifying the left-bottom-back corner and right-upper-front corner of a box with voxel size dx$\times$dy$\times$dz, (if only one value is given the voxels are cubic). We use a macroscopic volume average approach, where the volume fraction or porosity is the free space in which particles can diffuse (typically 0.2 in the brain) and the tortuosity is the increase in path length particles take due to obstacles.
k = rxd.Species(ecs, d=2.62, name='k', charge=1, initial=3)
The species diffusing in the ECS are define in the same was as intracellular species.
The extracellular concentrations are accessible by states3d
. Here we animate the concentrations during a 10s simulation.
from matplotlib import pyplot, animation
from IPython.display import HTML
# use a large timestep as this model only has slow diffusive dynamics.
h.dt = 10
#Create an animation of average (over depth) of the K+ concentration
def runsim(species, min_conc=3, max_conc=4, frames=1000):
h.finitialize()
fig = pyplot.figure()
im = pyplot.imshow(species[ecs].states3d.mean(2), vmin=min_conc, vmax=max_conc)
pyplot.axis('off')
def init():
im.set_data(species[ecs].states3d.mean(2))
return [im]
def animate(i):
h.fadvance()
im.set_data(species[ecs].states3d.mean(2))
return [im]
anim = animation.FuncAnimation(fig, animate, init_func=init, frames=frames, interval=10)
ret = HTML(anim.to_html5_video())
pyplot.close()
return ret
runsim(k, min_conc=3, max_conc=4)
Recording from a single voxel
We can also record from the extracellular concentration in the same way other NEURON states are recorded. Here we record from a extracellular node by access it by location or index into the states3d matrix.
# record from an extracellular nodes given by the x,y,z location
kecs_vec0 = h.Vector()
kecs_vec0.record(k[ecs].node_by_location(0, 0, 0)._ref_value) # same as k[ecs].node_by_ijk(20,15,7)
# record the same node by it's index into stated3d
kecs_vec1 = h.Vector()
kecs_vec1.record(k[ecs].node_by_ijk(22, 15, 7)._ref_value) # same as k[ecs].node_by_location(50,0,0)
# record the time
t_vec = h.Vector()
t_vec.record(h._ref_t)
# run the simulation
h.finitialize()
h.continuerun(1000)
# plot the concentations
pyplot.plot(t_vec, kecs_vec0, label='near the cell')
pyplot.plot(t_vec, kecs_vec1, label='far from the cell')
Boundary conditions
By default the extracellular region uses Neumann boundary conditions (where these is zero flux at the boundary), these means our simple example will eventually fill with K$^+$. One way to avoid this is to used Dirichlet boundary conditions where the concentration is set at the boundaries.
del k
k = rxd.Species(ecs, d=2.62, name='k', charge=1, initial=3, ecs_boundary_conditions=3)
runsim(k, min_conc=3, max_conc=4)
Extracellular reactions
Extracellular reactions are specified in the same way as intracellular ones. For example, this simple phenomenological model of potassium buffering by astrocytes.
from neuron.crxd import rxdmath
kb = 0.0008
kth = 15.0
kf = kb / (1.0 + rxdmath.exp(-(k - kth)/1.15))
Bmax = 10
A = rxd.Species(ecs, name='buffer', charge=1, d=0, initial=Bmax)
AK = rxd.Species(ecs, name='bound', charge=1, d=0, initial=0)
The buffering uses two additional species, A
is the unbound buffer and AK
is the bound buffer.
buffering = rxd.Reaction(k + A, AK, kf, kb)
The reaction is specified in the same way as intracellular rxd. Here we animate the concentration of bound buffer during a 10s simulation.
runsim(AK, min_conc=0, max_conc=0.1)
Inhomogeneities
Initial conditions
Initial condition in the ECS need not be a scalar concentration. Suppose we only want to apply the buffer shown above to a sphere in the middle of the ECS, this can be achieved by passing function that takes a node as an argument.
del A, buffering
A = rxd.Species(ecs,name='buffer', charge=1, d=0,
initial=lambda nd: 10 if (nd.x**2) + (nd.y**2) + nd.x**2 < 25**2 else 0)
buffering = rxd.Reaction(k + A, AK, kf, kb)
Anisotropy
The diffusion coefficient for a species can be different in each direction.
del A, AK, buffering, k
k = rxd.Species(ecs, d=(2.62, 2.62, 0), name='k', charge=1, initial=3)
Tissue characteristics
The volume faction and tortuosity of the extracellular space can also vary with location by passing a function coordinates x, y, z
. This is used in in ModelDB:238892 as a way to model edema following ischemic stroke.