-
Notifications
You must be signed in to change notification settings - Fork 1
/
HH2.mod
143 lines (121 loc) · 2.47 KB
/
HH2.mod
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
TITLE Hippocampal HH channels
:
: Fast Na+ and K+ currents responsible for action potentials
: Iterative equations
:
: Equations modified by Traub, for Hippocampal Pyramidal cells, in:
: Traub & Miles, Neuronal Networks of the Hippocampus, Cambridge, 1991
:
: range variable vtraub adjust threshold
:
: Written by Alain Destexhe, Salk Institute, Aug 1992
:
: Modified Oct 96 for compatibility with Windows: trap low values of arguments
:
INDEPENDENT {t FROM 0 TO 1 WITH 1 (ms)}
NEURON {
SUFFIX hh2
USEION na READ ena WRITE ina
USEION k READ ek WRITE ik
RANGE gnabar, gkbar, vtraub
RANGE m_inf, h_inf, n_inf
RANGE tau_m, tau_h, tau_n
RANGE m_exp, h_exp, n_exp
}
UNITS {
(mA) = (milliamp)
(mV) = (millivolt)
}
PARAMETER {
gnabar = .003 (mho/cm2)
gkbar = .005 (mho/cm2)
ena = 50 (mV)
ek = -90 (mV)
celsius = 36 (degC)
dt (ms)
v (mV)
vtraub = -63 (mV)
}
STATE {
m h n
}
ASSIGNED {
ina (mA/cm2)
ik (mA/cm2)
il (mA/cm2)
m_inf
h_inf
n_inf
tau_m
tau_h
tau_n
m_exp
h_exp
n_exp
tadj
}
BREAKPOINT {
SOLVE states
ina = gnabar * m*m*m*h * (v - ena)
ik = gkbar * n*n*n*n * (v - ek)
}
:DERIVATIVE states { : exact Hodgkin-Huxley equations
: evaluate_fct(v)
: m' = (m_inf - m) / tau_m
: h' = (h_inf - h) / tau_h
: n' = (n_inf - n) / tau_n
:}
PROCEDURE states() { : exact when v held constant
evaluate_fct(v)
m = m + m_exp * (m_inf - m)
h = h + h_exp * (h_inf - h)
n = n + n_exp * (n_inf - n)
VERBATIM
return 0;
ENDVERBATIM
}
UNITSOFF
INITIAL {
:
: Q10 was assumed to be 3 for both currents
:
tadj = 3.0 ^ ((celsius-36)/ 10 )
m = 0
h = 0
n = 0
}
PROCEDURE evaluate_fct(v(mV)) { LOCAL a,b,v2
v2 = v - vtraub : convert to traub convention
: a = 0.32 * (13-v2) / ( Exp((13-v2)/4) - 1)
a = 0.32 * vtrap(13-v2, 4)
: b = 0.28 * (v2-40) / ( Exp((v2-40)/5) - 1)
b = 0.28 * vtrap(v2-40, 5)
tau_m = 1 / (a + b) / tadj
m_inf = a / (a + b)
a = 0.128 * Exp((17-v2)/18)
b = 4 / ( 1 + Exp((40-v2)/5) )
tau_h = 1 / (a + b) / tadj
h_inf = a / (a + b)
: a = 0.032 * (15-v2) / ( Exp((15-v2)/5) - 1)
a = 0.032 * vtrap(15-v2, 5)
b = 0.5 * Exp((10-v2)/40)
tau_n = 1 / (a + b) / tadj
n_inf = a / (a + b)
m_exp = 1 - Exp(-dt/tau_m)
h_exp = 1 - Exp(-dt/tau_h)
n_exp = 1 - Exp(-dt/tau_n)
}
FUNCTION vtrap(x,y) {
if (fabs(x/y) < 1e-6) {
vtrap = y*(1 - x/y/2)
}else{
vtrap = x/(Exp(x/y)-1)
}
}
FUNCTION Exp(x) {
if (x < -100) {
Exp = 0
}else{
Exp = exp(x)
}
}