Parallel odor processing by mitral and middle tufted cells in the OB (Cavarretta et al 2016, 2018)

 Download zip file   Auto-launch 
Help downloading and running models
Accession:240116
"[...] experimental findings suggest that MC and mTC may encode parallel and complementary odor representations. We have analyzed the functional roles of these pathways by using a morphologically and physiologically realistic three-dimensional model to explore the MC and mTC microcircuits in the glomerular layer and deeper plexiform layers. [...]"
Reference:
1 . Cavarretta F, Burton SD, Igarashi KM, Shepherd GM, Hines ML, Migliore M (2018) Parallel odor processing by mitral and middle tufted cells in the olfactory bulb. Sci Rep 8:7625 [PubMed]
2 . Cavarretta F, Marasco A, Hines ML, Shepherd GM, Migliore M (2016) Glomerular and Mitral-Granule Cell Microcircuits Coordinate Temporal and Spatial Information Processing in the Olfactory Bulb. Front Comput Neurosci 10:67 [PubMed]
Citations  Citation Browser
Model Information (Click on a link to find other models with that property)
Model Type: Realistic Network;
Brain Region(s)/Organism: Olfactory bulb;
Cell Type(s): Olfactory bulb main tufted middle GLU cell; Olfactory bulb main interneuron granule MC GABA cell; Olfactory bulb main interneuron granule TC GABA cell; Olfactory bulb (accessory) mitral cell; Olfactory bulb main tufted cell external; Olfactory bulb short axon cell;
Channel(s): I A; I Na,t; I_Ks; I K;
Gap Junctions: Gap junctions;
Receptor(s): AMPA; GabaA; NMDA;
Gene(s):
Transmitter(s): Glutamate; Gaba;
Simulation Environment: NEURON;
Model Concept(s): Action Potentials; Action Potential Initiation; Active Dendrites; Long-term Synaptic Plasticity; Synaptic Integration; Synchronization; Pattern Recognition; Spatio-temporal Activity Patterns; Temporal Pattern Generation; Sensory coding; Sensory processing; Olfaction;
Implementer(s): Cavarretta, Francesco [francescocavarretta at hotmail.it]; Hines, Michael [Michael.Hines at Yale.edu];
Search NeuronDB for information about:  Olfactory bulb main interneuron granule MC GABA cell; Olfactory bulb main tufted middle GLU cell; Olfactory bulb main interneuron granule TC GABA cell; GabaA; AMPA; NMDA; I Na,t; I A; I K; I_Ks; Gaba; Glutamate;
/
modeldb-bulb3d
sim
ampanmda.mod
distrt.mod *
fi.mod
fi_stdp.mod *
gap.mod
Gfluct.mod
kamt.mod
kdrmt.mod
ks.mod
naxn.mod
orn.mod
ThreshDetect.mod *
all.py
all2all.py *
assembly.py
balance.py *
bindict.py
binsave.py
binspikes.py
blanes.hoc
blanes.py
blanes_exc_conn.txt
blanes6.dic
bulb3dtest.py
cancel.py
catfiles.sh
cellreader.py
cellwriter.py
cfg27.py
common.py
complexity.py *
convertdic.py
destroy_model.py
determine_connections.py
distribute.py *
dsac.py
Eta.txt *
fillgloms.py
fixnseg.hoc *
g_conn_stats.py
gapjunc.py
gen_weights.py
geodist.py
geodist.txt
getmitral.py
gidfunc.py
GJ.py
gj_nrn.hoc
Glom.py *
granule.hoc
granules.py
graphmeat.py
grow.py
growdef.py *
growout.py
job
Kod.txt *
lateral_connections.py
loadbalutil.py *
lpt.py *
mcgrow.py
MCrealSoma.py *
mgrs.py
misc.py
mitral.hoc
mkassembly.py
mkmitral.py
modeldata.py
mtgrow.py
MTrealSoma.py
MTrealSoma2.py
mtufted.hoc
multisplit_distrib.py
net_mitral_centric.py
Nod.txt *
odors.py
odorstim.py
odstim2.txt *
pad.txt *
params.py
parrun.py
pathdist.py
realgloms.txt *
runsim.py
spike2file.hoc *
spk2weight.py
split.py
subsetsim.py
test_complexity.py
txt2bin.py
util.py *
vrecord.py
weightsave.py
                            
TITLE simple NMDA receptors

: Hines combined AMPA and NMDA and spike dependent plasticity

: Modified from the original AMPA.mod, M.Migliore Jan 2003
: A weight of 0.0035 gives a peak conductance of 1nS in 0Mg

COMMENT
-----------------------------------------------------------------------------

	Simple model for glutamate AMPA receptors
	=========================================

  - FIRST-ORDER KINETICS, FIT TO WHOLE-CELL RECORDINGS

    Whole-cell recorded postsynaptic currents mediated by AMPA/Kainate
    receptors (Xiang et al., J. Neurophysiol. 71: 2552-2556, 1994) were used
    to estimate the parameters of the present model; the fit was performed
    using a simplex algorithm (see Destexhe et al., J. Computational Neurosci.
    1: 195-230, 1994).

  - SHORT PULSES OF TRANSMITTER (0.3 ms, 0.5 mM)

    The simplified model was obtained from a detailed synaptic model that 
    included the release of transmitter in adjacent terminals, its lateral 
    diffusion and uptake, and its binding on postsynaptic receptors (Destexhe
    and Sejnowski, 1995).  Short pulses of transmitter with first-order
    kinetics were found to be the best fast alternative to represent the more
    detailed models.

  - ANALYTIC EXPRESSION

    The first-order model can be solved analytically, leading to a very fast
    mechanism for simulating synapses, since no differential equation must be
    solved (see references below).



References

   Destexhe, A., Mainen, Z.F. and Sejnowski, T.J.  An efficient method for
   computing synaptic conductances based on a kinetic model of receptor binding
   Neural Computation 6: 10-14, 1994.  

   Destexhe, A., Mainen, Z.F. and Sejnowski, T.J. Synthesis of models for
   excitable membranes, synaptic transmission and neuromodulation using a 
   common kinetic formalism, Journal of Computational Neuroscience 1: 
   195-230, 1994.


-----------------------------------------------------------------------------
ENDCOMMENT



NEURON {
	POINT_PROCESS AmpaNmda
	RANGE R, g, mg, inmda, iampa, gnmda, gampa
	RANGE x, mgid, ggid, srcgid, gmax, training
	NONSPECIFIC_CURRENT i
	GLOBAL Cdur, Alpha, Beta, E, Rinf, Rtau, ampatau
	GLOBAL gampafactor, nmdafactor
	GLOBAL ltdinvl, ltpinvl, sighalf, sigslope, sigexp
}
UNITS {
	(nA) = (nanoamp)
	(mV) = (millivolt)
	(umho) = (micromho)
	(mM) = (milli/liter)
}

PARAMETER {

	Cdur	= 1		(ms)	: transmitter duration (rising phase)
	Alpha	= 0.35		(/ms)	: forward (binding) rate
	Beta	= 0.012		(/ms)	: backward (unbinding) rate
	E	= 0	(mV)		: reversal potential
	mg	= 1    (mM)		: external magnesium concentration
	gmax = 2 (umho)		: normally 2
	gampafactor = 0.001 (1)
	nmdafactor = 0.07 (1)
	ltdinvl = 250 (ms)		: longer intervals, no change
	ltpinvl = 33.33 (ms)		: shorter interval, LTP
	sighalf = 25 (1)
	sigslope = 5 (1)
        sigexp = 4
	ampatau = 3 (ms)
	x = 0 (um) : cartesian synapse location
	mgid = -1 : associated mitral gid
	ggid = -1 : associated granule gid
	srcgid = -1 : gid of the mitral detector
	training = 1 : is on
}


ASSIGNED {
	v		(mV)		: postsynaptic voltage
	i 		(nA)		: total current = iampa+inmda
	inmda 		(nA)		: current = gnmda*(v - E)
	iampa 		(nA)		: current = gampa*(v - E)
	gnmda 		(umho)		: 
	Rinf				: steady state channels open
	Rtau		(ms)		: time constant of channel binding
	synon
}

STATE {Ron Roff
	gampa 		(umho)
}

INITIAL {
	PROTECT Rinf = Alpha / (Alpha + Beta)
	PROTECT Rtau = 1 / (Alpha + Beta)
	synon = 0
	gampa = 0
}

BREAKPOINT {
	SOLVE release METHOD cnexp
	gnmda = mgblock(v)*(Ron + Roff)*gmax*nmdafactor
	inmda = gnmda*(v - E)
	iampa = gampa*(v - E)
	i = iampa + inmda
}

DERIVATIVE release {
	Ron' = (synon*Rinf - Ron)/Rtau
	Roff' = -Beta*Roff
	gampa' = -gampa/ampatau
}

: following supports both saturation from single input and
: summation from multiple inputs
: if spike occurs during CDur then new off time is t + CDur
: ie. transmitter concatenates but does not summate
: Note: automatic initialization of all reference args to 0 except first


FUNCTION mgblock(v(mV)) {
	TABLE 
	DEPEND mg
	FROM -140 TO 80 WITH 1000

	: from Jahr & Stevens

	mgblock = 1 / (1 + exp(0.062 (/mV) * -v) * (mg / 3.57 (mM)))
}

FUNCTION plast(step(1))(1) {
	plast = (1 - 1/(1 + exp((step - sighalf)/sigslope)))^sigexp

}

NET_RECEIVE(weight, s, w, tlast (ms), r0, t0 (ms)) {
	INITIAL {
		:s = 0 
                if(s == 0) {
                  w = 0
                } else {
		  w = weight*plast(s)
                }
		tlast = -1e9 (ms)
		r0 = 0
		t0 = -1e9 (ms)
	}
	: flag is an implicit argument of NET_RECEIVE and  normally 0
        if (flag == 0) { : a spike, so turn on if not already in a Cdur pulse
		: plasticity affects this spike. If desired to affect
		: the next spike then put following group after
		: net_send
		if (t - tlast < ltpinvl) { : LTP
                        if (training) {
			  s = s + 1
			  :if (s > 2*sighalf) { s = 2*sighalf }
                          if(s>75) {s=75}
                        }
		}else if (t - tlast > ltdinvl) { : no change
		}else{ : LTD
                      if (training) {
			s = s - 1
			if (s < 0) { s = 0 }
                      }
		}
                         
		tlast = t
                         
                if(s == 0) {
                    w = 0
                } else {
                    w = weight*plast(s)
                }
                         
		gampa = gampa + w*gmax*gampafactor
		r0 = r0*exp(-Beta*(t - t0))
		t0 = t
		synon = synon + w
		Ron = Ron + r0
		Roff = Roff - r0
		: come again in Cdur with flag = current value of w+1
		net_send(Cdur, w + 1)
        }else{ : turn off what was added Cdur ago
		r0 = (flag-1)*Rinf + (r0 - (flag-1)*Rinf)*exp(-(t - t0)/Rtau)
		t0 = t
		synon = synon - (flag-1)
		Ron = Ron - r0
		Roff = Roff + r0
	}
}