Computational analysis of NN activity and spatial reach of sharp wave-ripples (Canakci et al 2017)

 Download zip file 
Help downloading and running models
Accession:230861
Network oscillations of different frequencies, durations and amplitudes are hypothesized to coordinate information processing and transfer across brain areas. Among these oscillations, hippocampal sharp wave-ripple complexes (SPW-Rs) are one of the most prominent. SPW-Rs occurring in the hippocampus are suggested to play essential roles in memory consolidation as well as information transfer to the neocortex. To-date, most of the knowledge about SPW-Rs comes from experimental studies averaging responses from neuronal populations monitored by conventional microelectrodes. In this work, we investigate spatiotemporal characteristics of SPW-Rs and how microelectrode size and distance influence SPW-R recordings using a biophysical model of hippocampus. We also explore contributions from neuronal spikes and synaptic potentials to SPW-Rs based on two different types of network activity. Our study suggests that neuronal spikes from pyramidal cells contribute significantly to ripples while high amplitude sharp waves mainly arise from synaptic activity. Our simulations on spatial reach of SPW-Rs show that the amplitudes of sharp waves and ripples exhibit a steep decrease with distance from the network and this effect is more prominent for smaller area electrodes. Furthermore, the amplitude of the signal decreases strongly with increasing electrode surface area as a result of averaging. The relative decrease is more pronounced when the recording electrode is closer to the source of the activity. Through simulations of field potentials across a high-density microelectrode array, we demonstrate the importance of finding the ideal spatial resolution for capturing SPW-Rs with great sensitivity. Our work provides insights on contributions from spikes and synaptic potentials to SPW-Rs and describes the effect of measurement configuration on LFPs to guide experimental studies towards improved SPW-R recordings.
Reference:
1 . Canakci S, Toy MF, Inci AF, Liu X, Kuzum D (2017) Computational analysis of network activity and spatial reach of sharp wave-ripples. PLoS One 12:e0184542 [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: Hippocampus;
Cell Type(s): Hippocampus CA1 pyramidal GLU cell; Hippocampus CA1 basket cell;
Channel(s): I Na,t; I A; I K; I h;
Gap Junctions: Gap junctions;
Receptor(s): NMDA; GabaA; Glutamate;
Gene(s):
Transmitter(s):
Simulation Environment: NEURON;
Model Concept(s): Oscillations; Spatio-temporal Activity Patterns;
Implementer(s): Canakci, Sadullah [scanakci at bu.edu]; Inci, Ahmet F [afinci at sabanciuniv,edu]; Toy, Faruk [faruk.toy at metu.edu.tr]; Liu, Xin [xil432 at end.ucsd.edu]; Kuzum, Duygu [dkuzum at eng.ucsd.edu];
Search NeuronDB for information about:  Hippocampus CA1 pyramidal GLU cell; GabaA; NMDA; Glutamate; I Na,t; I A; I K; I h;
/
CanakciEtAl2017
mod
ca1ih.mod *
ca1ika.mod *
ca1ikdr.mod *
ca1ina.mod *
caolmw.mod *
capr.mod *
edrive.mod *
exp2synampa.mod *
exp2synampapre.mod *
exp2synnmda.mod *
exp2synnmdaperm.mod *
exp2synnmdapre.mod *
fvpre.mod *
gap.mod *
halfgap.mod *
icaolmw.mod *
icapr.mod *
iholmkop.mod *
iholmw.mod *
ihpyrkop.mod *
kahppr.mod *
kaolmkop.mod *
kapyrkop.mod *
kcaolmw.mod *
kcpr.mod *
kdrbwb.mod *
kdrca1.mod *
kdrolmkop.mod *
kdrpr.mod *
kdrpyrkop.mod *
na3n.mod *
nafbwb.mod *
nafolmkop.mod *
nafpr.mod *
nafpyrkop.mod *
noisesyn.mod
par_ggap.mod *
vecevent.mod *
xtra.mod *
aux_fun.inc *
                            
TITLE na3
: Na current 
: modified from Jeff Magee. M.Migliore may97
: added sh to account for higher threshold M.Migliore, Apr.2002

NEURON {
	SUFFIX na3
	USEION na READ ena WRITE ina
	RANGE  gbar, ar, sh
	GLOBAL minf, hinf, mtau, htau, sinf, taus,qinf, thinf
}

PARAMETER {
	sh   = 8	(mV)
	gbar = 0.032   	(mho/cm2)	
								
	tha  =  -30	(mV)		: v 1/2 for act	
	qa   = 7.2	(mV)		: act slope (4.5)		
	Ra   = 0.4	(/ms)		: open (v)		
	Rb   = 0.124 	(/ms)		: close (v)		

	thi1  = -45	(mV)		: v 1/2 for inact 	
	thi2  = -45 	(mV)		: v 1/2 for inact 	
	qd   = 1.5	(mV)	        : inact tau slope
	qg   = 1.5      (mV)
	mmin=0.02	
	hmin=0.5			
	q10=2
	Rg   = 0.01 	(/ms)		: inact recov (v) 	
	Rd   = .03 	(/ms)		: inact (v)	
	qq   = 10        (mV)
	tq   = -55      (mV)

	thinf  = -50 	(mV)		: inact inf slope	
	qinf  = 4 	(mV)		: inact inf slope 

        vhalfs=-60	(mV)		: slow inact.
        a0s=0.0003	(ms)		: a0s=b0s
        zetas=12	(1)
        gms=0.2		(1)
        smax=10		(ms)
        vvh=-58		(mV) 
        vvs=2		(mV)
        ar=0.5		(1)		: 1=no inact., 0=max inact.
	ena = 55		(mV)            : must be explicitly def. in hoc
	celsius
	v 		(mV)
}


UNITS {
	(mA) = (milliamp)
	(mV) = (millivolt)
	(pS) = (picosiemens)
	(um) = (micron)
} 

ASSIGNED {
	ina 		(mA/cm2)
	thegna		(mho/cm2)
	minf 		hinf 		
	mtau (ms)	htau (ms) 	
	sinf (ms)	taus (ms)
}
 

STATE { m h s}

BREAKPOINT {
        SOLVE states METHOD cnexp
        thegna = gbar*m*m*m*h*s
	ina = thegna * (v - ena)
} 

INITIAL {
	trates(v,ar,sh)
	m=minf  
	h=hinf
	s=sinf
}


FUNCTION alpv(v(mV)) {
         alpv = 1/(1+exp((v-vvh-sh)/vvs))
}
        
FUNCTION alps(v(mV)) {  
  alps = exp(1.e-3*zetas*(v-vhalfs-sh)*9.648e4/(8.315*(273.16+celsius)))
}

FUNCTION bets(v(mV)) {
  bets = exp(1.e-3*zetas*gms*(v-vhalfs-sh)*9.648e4/(8.315*(273.16+celsius)))
}

LOCAL mexp, hexp, sexp

DERIVATIVE states {   
        trates(v,ar,sh)      
        m' = (minf-m)/mtau
        h' = (hinf-h)/htau
        s' = (sinf - s)/taus
}

PROCEDURE trates(vm,a2,sh2) {  
        LOCAL  a, b, c, qt
        qt=q10^((celsius-24)/10)
	a = trap0(vm,tha+sh2,Ra,qa)
	b = trap0(-vm,-tha-sh2,Rb,qa)
	mtau = 1/(a+b)/qt
        if (mtau<mmin) {mtau=mmin}
	minf = a/(a+b)

	a = trap0(vm,thi1+sh2,Rd,qd)
	b = trap0(-vm,-thi2-sh2,Rg,qg)
	htau =  1/(a+b)/qt
        if (htau<hmin) {htau=hmin}
	hinf = 1/(1+exp((vm-thinf-sh2)/qinf))
	c=alpv(vm)
        sinf = c+a2*(1-c)
        taus = bets(vm)/(a0s*(1+alps(vm)))
        if (taus<smax) {taus=smax}
}

FUNCTION trap0(v,th,a,q) {
	if (fabs(v-th) > 1e-6) {
	        trap0 = a * (v - th) / (1 - exp(-(v - th)/q))
	} else {
	        trap0 = a * q
 	}
}