ModelDB is moving. Check out our new site at https://modeldb.science. The corresponding page is https://modeldb.science/225080.

CA1 pyr cell: Inhibitory modulation of spatial selectivity+phase precession (Grienberger et al 2017)

 Download zip file 
Help downloading and running models
Accession:225080
Spatially uniform synaptic inhibition enhances spatial selectivity and temporal coding in CA1 place cells by suppressing broad out-of-field excitation.
Reference:
1 . Grienberger C, Milstein AD, Bittner KC, Romani S, Magee JC (2017) Inhibitory suppression of heterogeneously tuned excitation enhances spatial coding in CA1 place cells. Nat Neurosci 20:417-426 [PubMed]
Model Information (Click on a link to find other models with that property)
Model Type: Neuron or other electrically excitable cell; Realistic Network;
Brain Region(s)/Organism: Hippocampus;
Cell Type(s): Hippocampus CA1 pyramidal GLU cell;
Channel(s):
Gap Junctions:
Receptor(s): NMDA;
Gene(s):
Transmitter(s):
Simulation Environment: NEURON; Python;
Model Concept(s): Active Dendrites; Detailed Neuronal Models; Place cell/field; Synaptic Integration; Short-term Synaptic Plasticity; Spatial Navigation; Feature selectivity;
Implementer(s): Milstein, Aaron D. [aaronmil at stanford.edu];
Search NeuronDB for information about:  Hippocampus CA1 pyramidal GLU cell; NMDA;
/
GrienbergerEtAl2017
morphologies
readme.txt
ampa_kin.mod *
exp2EPSC.mod
exp2EPSG.mod
exp2EPSG_NMDA.mod
gaba_a_kin.mod *
h.mod
kad.mod *
kap.mod *
kdr.mod *
km2.mod
nas.mod
nax.mod
nmda_kin2.mod
nmda_kin3.mod
nmda_kin5.mod *
pr.mod *
vecevent.mod *
batch_EPSP_attenuation.sh
batch_place_cell_r_inp.sh
batch_place_cell_record_i_syn.sh
batch_place_cell_single_compartment.sh
batch_place_cell_subtr_inh.sh
batch_place_cell_subtr_inh_shifted.sh
batch_place_cell_subtr_inh_vclamp.sh
batch_process_i_syn_files.sh
batch_rinp.sh
batch_spine_attenuation_ratio.sh
build_expected_EPSP_reference.sh
build_expected_EPSP_reference_controller.py
build_expected_EPSP_reference_engine.py
consolidate_i_syn_files.py
consolidate_tracked_spine_data.py
fit_parameter_exponential_distribution.py
function_lib.py
optimize_AMPA_KIN.py
optimize_dendritic_excitability_020416.py
optimize_GABA_A_KIN.py
optimize_NMDA_KIN2.py
parallel_branch_cooperativity.sh
parallel_branch_cooperativity_no_nmda.sh
parallel_clustered_branch_cooperativity_nmda_controller_110315.py
parallel_clustered_branch_cooperativity_nmda_engine_110315.py
parallel_EPSP_attenuation_controller.py
parallel_EPSP_attenuation_engine.py
parallel_EPSP_i_attenuation_controller.py
parallel_EPSP_i_attenuation_engine.py
parallel_expected_EPSP_controller.py
parallel_expected_EPSP_engine.py
parallel_optimize_branch_cooperativity.sh
parallel_optimize_branch_cooperativity_nmda_kin3_controller.py
parallel_optimize_branch_cooperativity_nmda_kin3_engine.py
parallel_optimize_EPSP_amp_controller.py
parallel_optimize_EPSP_amp_engine.py
parallel_optimize_pr.sh
parallel_optimize_pr_controller_020116.py
parallel_optimize_pr_engine_020116.py
parallel_rinp_controller.py
parallel_rinp_engine.py
parallel_spine_attenuation_ratio_controller.py
parallel_spine_attenuation_ratio_engine.py
plot_channel_distributions.py
plot_NMDAR_kinetics.py
plot_results.py
plot_spine_traces.py
plot_synaptic_conductance_facilitation.py
process_i_syn_files.py
record_bAP_attenuation.py
simulate_place_cell_no_precession.py
simulate_place_cell_single_compartment.py
simulate_place_cell_single_compartment_no_nmda.py
simulate_place_cell_subtr_inh.py
simulate_place_cell_subtr_inh_add_noise.py
simulate_place_cell_subtr_inh_add_noise_no_na.py
simulate_place_cell_subtr_inh_no_na.py
simulate_place_cell_subtr_inh_no_nmda_no_na.py
simulate_place_cell_subtr_inh_r_inp.py
simulate_place_cell_subtr_inh_rec_i_syn.py
simulate_place_cell_subtr_inh_shifted.py
simulate_place_cell_subtr_inh_silent.py
simulate_place_cell_subtr_inh_vclamp.py
specify_cells.py
                            
COMMENT
-----------------------------------------------------------------------------

GluN-R (NMDA-type Glutamate Receptors)

Postsynaptic current is additionally constrained by voltage-dependent channel
block by extracellular Mg.

-----------------------------------------------------------------------------
Synaptic mechanism based on a simplified model of transmitter binding to
postsynaptic receptors.

Modified by Aaron Milstein, 2015.

Modification of original code by:

A. Destexhe & Z. Mainen, The Salk Institute, March 12, 1993.
Last modif. Sept 8, 1993.

Reference:

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

-----------------------------------------------------------------------------

Upon arrival of a presynaptic spike (a net_event), the concentration of the
neurotransmitter C in the synaptic cleft is briefly stepped to concentration
Cmax for duration Cdur.

C     _____ . . . . . . Cmax
      |   |
 _____|   |______ . . . 0
     t0   t0 + Cdur

The receptors then bind transmitter, change conformation, and open their
channel according to the following kinetic scheme:

       ---[C] * kon-->        -------CC----->        -----Beta----->
C + Ru <-----koff-----   Rb   <------CO------   Rc   <----Alpha-----   Ro

where Ru, Rb, Rc, and Ro are respectively the fraction of channels in the
unbound, closed bound, closed cleft, and open states of the postsynaptic
receptor. kon and koff are the binding and unbinding rate constants, CC and
CO are the closing and opening rates of the receptor ligand-binding cleft,
and Beta and Alpha are the opening and closing rates of the channel.

The maximal conductance of the channel can be set for each instance of the
synaptic mechanism by specifying gmax, and the relative weight of events
can be set through the weight parameter of the associated netcon object.

The postsynaptic current is given by:

i = weight * gmax * Ro * (V-Erev)

-----------------------------------------------------------------------------

ENDCOMMENT

NEURON {
	POINT_PROCESS NMDA_KIN3
	RANGE Cmax, Cdur, kon, koff, CC, CO, Beta, Alpha, Erev, Kd, gamma, mg, gmax, g, B, kin_scale
	NONSPECIFIC_CURRENT i
}
UNITS {
	(nA) = (nanoamp)
	(mV) = (millivolt)
	(umho) = (micromho)
	(mM) = (milli/liter)
}

PARAMETER {

	Cmax        = 1.        (mM)        : transmitter concentration during release event
	Cdur	    = 0.3       (ms)		: transmitter duration (rising phase)
	kon         = 86.89     (/mM/ms)    : unbound receptor ligand-binding rate
    koff        = 0.69      (/ms)       : bound receptor ligand-unbinding rate
    CC          = 9.64      (/ms)       : bound receptor cleft closing rate
    CO          = 2.60      (/ms)       : bound receptor cleft opening rate
    Beta	    = 0.68      (/ms)	    : channel opening rate
    Alpha       = 0.079     (/ms)       : open channel closing rate
	Erev	    = 0.        (mV)		: reversal potential
	Kd          = 7.51      (mM)        : modulate Mg concentration dependence
    gamma       = 0.100     (/mV)       : modulate slope of Mg sensitivity
    mg          = 1.0       (mM)        : extracellular Mg concentration
    gmax	    = 0.00361   (umho)	    : maximum conductance
    kin_scale   = 1.81      (1)         : scale voltage sensitivity of decay kinetics
}


ASSIGNED {
	v		(mV)		: postsynaptic voltage
	i 		(nA)		: current = g * (v - Erev)
	g 		(umho)		: conductance
	C       (mM)        : unbound transmitter concentration
    B                   : fraction of channels not blocked by extracellular Mg
    kin_factor          : exponential scaling of decay rates with voltage
    scale               : allow netcon weight to scale conductance
}

STATE {
    Ru                  : fraction of receptors not bound to transmitter
    Rb                  : fraction of receptors bound to transmitter
    Rc                  : fraction of receptors in closed cleft state
    Ro                  : fraction of channels in open state
}

INITIAL {
	C = 0.
    Ru = 1.
    Rb = 0.
    Rc = 0.
    Ro = 0.
    scale = 1.
    mgblock(v)
}

BREAKPOINT {
    mgblock(v)
    SOLVE kstates METHOD sparse
	g = scale * gmax * B * Ro
    i = g * (v - Erev)
}

KINETIC kstates {
    ~ Ru <-> Rb     (C * kon, koff)
    ~ Rb <-> Rc     (CC, CO * kin_factor)
    ~ Rc <-> Ro     (Beta, Alpha * kin_factor)
}

NET_RECEIVE(weight) {
    if (flag == 0) { : a new spike received
        C = Cmax
        scale = weight
        net_send(Cdur, 1)
    } else {    : a self event
        C = 0.
    }
}

PROCEDURE mgblock(v(mV)) {
	: from Jahr & Stevens
    B = 1. / (1. + exp(gamma * (-v)) * (mg / Kd))
    kin_factor = (1 - kin_scale) * B + kin_scale
}

Loading data, please wait...