Library

Parameters

SynapseElife.PreSynapseParamsType
struct PreSynapseParams

Presynaptic parameters

  • Firing events are processed separately from the main simulation (at src/OnlyStp.jl) it takes the firing structure as input from the function firingPattern.
  • Using the presynaptic stimulation times, the vesicle depletion and AP induced by EPSP are estimated, however one can use a tag to deactivate it (covering sub-threshold EPSP cases) as in dataProtocol.
  • The presynaptic part of our model is phenomenological, for instance, the variable Soma in src/OnlyStp.jl:38 was made to represent the voltage and can accumulate for faster frequencies but has an abstract unit.

Equations

based on D. Sterratt et al book; Principles of Computational Modelling in Neuroscience page 183

raterec = (R_0 - R) ⋅ τ_R ⋅ rrp

raterrp = (D_0 - D) ⋅ τ_D ⋅ rec

rateref = (R_0 - R) ⋅ ref_dt

Fields

  • τ_rec: recovery constant of pre calcium decay function Default: 20000

  • δ_ca: fraction of decay constant of pre calcium decay f Default: 0.0004

  • τ_pre: decay time constant of pre calcium Default: 20

  • τ_V: decay time constant for AP induced by EPSP Default: 40

  • δ_delay_AP: delay to EPSPs onset and evoked AP Default: 15.0

  • D_0: initial conditions ready releaseble pool Default: 25

  • R_0: initial conditions recovery pool Default: 30

  • τ_R: rate for D -> R Default: 5000

  • τ_D: rate for R -> D Default: 45000

  • τ_R_ref: rate for infinite reservoir -> R Default: 40000

  • s: sigmoid parameter for release probability Default: 2.0

  • h: sigmoid parameter for release probability Default: 0.7

  • sampling_rate: sampling rate for plotting / printing Default: 1.0

source
SynapseElife.SynapseParamsType
struct SynapseParams{Tp}

Postsynaptic parameters.


Units

  • time: ms
  • length: µm (area µm^2, volume µm^3)
  • voltage: mV
  • current: pA
  • conductance: nS
  • resistance: GOhm
  • capacitance: pF ( ms.pA/mV = ms.nS = ms/GOhm)
  • concentration: µM

Fields

  • LTP_region: Region LTP Default: VPolygon([[6.35, 1.4], [10, 1.4], [6.35, 29.5], [10, 29.5]])

  • LTD_region: Region LTD Default: VPolygon([[6.35, 1.4], [6.35, 23.25], [6.35, 29.5], [1.85, 11.327205882352938], [1.85, 23.25], [3.7650354609929075, 1.4], [5.650675675675676, 29.5]])

  • a_D: ACTIVATION RATES FOR THRESHOLDS Default: 0.1

  • b_D: Default: 2.0e-5

  • a_P: Default: 0.2

  • b_P: Default: 0.0001

  • t_D: Default: 18000

  • t_P: Default: 13000

  • K_D: Sigmoids controlling the rate of plasticity change Default: 80000.0

  • K_P: Default: 13000.0

  • rest_plstcty: Plasticity states Default: 100

  • t_end: SIMULATION Default: 100

  • sampling_rate: Default: 10

  • temp_rates: BIOPHYSICAL AND GHK PARAMETERS Default: 35.0

  • age: Default: 60.0

  • faraday: Default: 0.096485 * 0.001

  • Ca_ext: Default: 2500.0

  • Ca_infty: Default: 0.05

  • tau_ca: Default: 10.0

  • D_Ca: Default: 0.3338

  • f_Ca: Default: 0.1

  • perm: Default: -0.04583333333333333

  • z: Default: 2.0

  • gas: Default: 8.314e-6

  • p_release: Default: (0.004225803293622208, 1708.4124496514878, 1.3499793762587964, 0.6540248201173222)

  • trec: BACKPROPAGATION ATTENUATION Default: 2000

  • trec_soma: Default: 500

  • delta_decay: Default: 1.7279e-5

  • p_age_decay_bap: Default: (0.13525468256031167, 16.482800452454164, 5.564691354645679)

  • delta_soma: Default: 2.5e-5 * (pagedecaybap[3] / (1 + exp(pagedecaybap[1] * (age - pagedecay_bap[2]))))

  • delta_aux: Default: 2.304e-5

  • injbap: Default: 2.0

  • soma_dist: Default: 200.0

  • p_dist: Default: (0.019719018173341547, 230.3206470553394, 1.4313810030893268, 0.10406540965358434)

  • ϕ_dist: Default: pdist[4] + pdist[3] / (1 + exp(pdist[1] * (somadist - p_dist[2])))

  • I_clamp: Default: 0.0

  • gamma_Na: Na and K Default: 800.0

  • Erev_Na: Default: 50.0

  • gamma_K: Default: 40.0

  • Erev_K: Default: -90.0

  • p_nmda_frwd: NMDAr temperature modification Default: (-0.09991802053299291, -37.63132907014948, 1239.0673283348326, -1230.6805720050966)

  • frwd_T_chng_NMDA: Default: pnmdafrwd[4] + pnmdafrwd[3] / (1 + exp(pnmdafrwd[1] * (temprates - pnmda_frwd[2])))

  • p_nmda_bcwd: Default: (-0.10605060141396823, 98.99939433046647, 1621.6168608608068, 3.0368551011554143)

  • bcwd_T_chng_NMDA: Default: pnmdabcwd[4] + pnmdabcwd[3] / (1 + exp(pnmdabcwd[1] * (temprates - pnmda_bcwd[2])))

  • NMDA_N2A_ka: NMDAr kinetics (GluN2A type) Default: frwdTchng_NMDA * 34.0 * 0.001

  • NMDA_N2A_kb: Default: frwdTchng_NMDA * 17.0 * 0.001

  • NMDA_N2A_kc: Default: frwdTchng_NMDA * 127.0 * 0.001

  • NMDA_N2A_kd: Default: frwdTchng_NMDA * 580.0 * 0.001

  • NMDA_N2A_ke: Default: frwdTchng_NMDA * 2508.0 * 0.001

  • NMDA_N2A_kf: Default: frwdTchng_NMDA * 3449.0 * 0.001

  • NMDA_N2A_k_f: Default: bcwdTchng_NMDA * 662.0 * 0.001

  • NMDA_N2A_k_e: Default: bcwdTchng_NMDA * 2167.0 * 0.001

  • NMDA_N2A_k_d: Default: bcwdTchng_NMDA * 2610.0 * 0.001

  • NMDA_N2A_k_c: Default: bcwdTchng_NMDA * 161.0 * 0.001

  • NMDA_N2A_k_b: Default: bcwdTchng_NMDA * 120.0 * 0.001

  • NMDA_N2A_k_a: Default: bcwdTchng_NMDA * 60.0 * 0.001

  • NMDA_N2B_sa: NMDAr kinetics (GluN2B type) Default: frwdTchng_NMDA * 0.25 * 34.0 * 0.001

  • NMDA_N2B_sb: Default: frwdTchng_NMDA * 0.25 * 17.0 * 0.001

  • NMDA_N2B_sc: Default: frwdTchng_NMDA * 0.25 * 127.0 * 0.001

  • NMDA_N2B_sd: Default: frwdTchng_NMDA * 0.25 * 580.0 * 0.001

  • NMDA_N2B_se: Default: frwdTchng_NMDA * 0.25 * 2508.0 * 0.001

  • NMDA_N2B_sf: Default: frwdTchng_NMDA * 0.25 * 3449.0 * 0.001

  • NMDA_N2B_s_f: Default: bcwdTchng_NMDA * 0.23 * 662.0 * 0.001

  • NMDA_N2B_s_e: Default: bcwdTchng_NMDA * 0.23 * 2167.0 * 0.001

  • NMDA_N2B_s_d: Default: bcwdTchng_NMDA * 0.23 * 2610.0 * 0.001

  • NMDA_N2B_s_c: Default: bcwdTchng_NMDA * 0.23 * 161.0 * 0.001

  • NMDA_N2B_s_b: Default: bcwdTchng_NMDA * 0.23 * 120.0 * 0.001

  • NMDA_N2B_s_a: Default: bcwdTchng_NMDA * 0.23 * 60.0 * 0.001

  • p_nmda: Default: (0.004477162852447629, 2701.3929349701334, 58.38819453272428, 33.949463268365555)

  • gamma_nmda: Default: (pnmda[4] + pnmda[3] / (1 + exp(pnmda[1] * (Caext - p_nmda[2])))) * 0.001

  • p_age: Default: (0.09993657672916968, 25.102347872464193, 0.9642137892004939, 0.5075183905839776)

  • r_NMDA_age: Ratio N2B/N2A Default: rand(Normal(0, 0.05)) + page[4] + page[3] / (1 + exp(page[1] * (age - page[2])))

  • N_NMDA: Default: 15

  • N_N2B: Default: round((NNMDA * rNMDAage) / (rNMDA_age + 1))

  • N_N2A: Default: round(NNMDA / (rNMDA_age + 1))

  • Erev_nmda: Other NMDAr parameters Default: 0.0

  • Mg: Default: 1.0

  • p_ampa_frwd: AMPAr temperature modification Default: (-0.4737773089201679, 31.7248285571622, 10.273135485873242)

  • frwd_T_chng_AMPA: Default: pampafrwd[3] / (1 + exp(pampafrwd[1] * (temprates - pampa_frwd[2])))

  • p_ampa_bcwd: Default: (-0.36705555170278986, 28.976662403966674, 5.134547217640794)

  • bcwd_T_chng_AMPA: Default: pampabcwd[3] / (1 + exp(pampabcwd[1] * (temprates - pampa_bcwd[2])))

  • AMPA_k1: ** AMPAr kinetics Default: frwdTchng_AMPA * 1.6 * 1.0e7 * 1.0e-6 * 0.001

  • AMPA_k_1: Default: bcwdTchng_AMPA * 7400 * 0.001

  • AMPA_k_2: Default: bcwdTchng_AMPA * 0.41 * 0.001

  • AMPA_alpha: Default: 2600 * 0.001

  • AMPA_beta: Default: 9600 * 0.001

  • AMPA_delta_1: Default: 1500 * 0.001

  • AMPA_gamma_1: Default: 9.1 * 0.001

  • AMPA_delta_2: Default: 170 * 0.001

  • AMPA_gamma_2: Default: 42 * 0.001

  • AMPA_delta_0: Default: 0.003 * 0.001

  • AMPA_gamma_0: Default: 0.83 * 0.001

  • gamma_ampa1: AMPAr conductances Default: 0.5 * 0.031

  • gamma_ampa2: Default: 0.5 * 0.052

  • gamma_ampa3: Default: 0.5 * 0.073

  • N_ampa: Default: 120

  • Erev_ampa: Default: 0.0

  • N_GABA: GABAr Default: 34

  • p_Cl: Default: (0.09151696057098718, 0.6919298240788684, 243.5159017060495, -92.6496083089155)

  • Erev_Cl: GABAr, Cl reversal potential Default: pCl[4] + pCl[3] / (1 + exp(pCl[1] * (age - pCl[2])))

  • gamma_GABA: Default: 0.035

  • GABA_r_b1: Default: 1.0e6 * 1.0e-6 * 0.001 * 20

  • GABA_r_u1: Default: 1000.0 * 0.0046

  • GABA_r_b2: Default: 1.0e6 * 1.0e-6 * 0.001 * 10

  • GABA_r_u2: Default: 1000.0 * 0.0092

  • GABA_r_ro1: Default: 1000.0 * 0.0033

  • GABA_r_ro2: Default: 1000.0 * 0.0106

  • p_GABA: Default: (0.19127068198185954, 32.16771140618756, -1.2798050197287802, 1.470692263981145)

  • GABA_r_c1: Default: (pGABA[4] + pGABA[3] / (1 + exp(pGABA[1] * (temprates - p_GABA[2])))) * 1000.0 * 0.0098

  • GABA_r_c2: Default: (pGABA[4] + pGABA[3] / (1 + exp(pGABA[1] * (temprates - p_GABA[2])))) * 0.4

  • E_leak: Passive electrical properties Default: -70.0

  • g_leak: Default: 4.0e-6

  • Cm: Default: 0.006

  • R_a: Default: 0.01

  • D_dend: ** MORPHOLOGY: Dendritic properties** Default: 2.0

  • L_dend: Default: 1400

  • A_dend: Default: 2 * pi * (Ddend / 2) * Ldend

  • Vol_dend: Default: pi * (Ddend / 2) ^ 2 * Ldend

  • Cdend: Default: Cm * A_dend

  • CS_dend: Default: pi * (D_dend / 2) .^ 2

  • g_leakdend: Default: gleak * Adend

  • D_soma: MORPHOLOGY: Soma properties Default: 30

  • A_soma: Default: pi * D_soma ^ 2

  • Csoma: Default: Cm * A_soma

  • CS_soma: Default: pi * (D_soma / 2) .^ 2

  • g_leaksoma: Default: 15.0

  • g_diff: Default: Ddend / (4Ra)

  • Vol_sp: Spine properties Default: 0.03

  • A_sp: Default: 4 * pi * ((3Vol_sp) / (4pi)) ^ (2.0 / 3.0)

  • Csp: Default: Cm * A_sp

  • g_leaksp: Default: gleak * Asp

  • D_neck: Neck properties Default: 0.1

  • L_neck: Default: 0.2

  • CS_neck: Default: pi * (D_neck / 2) .^ 2

  • g_neck: Default: CSneck / (Lneck * R_a)

  • tau_diff: Default: Volsp / (2 * DCa * Dneck) + Lneck ^ 2 / (2D_Ca)

  • glu_width: SYNAPTIC GLUTAMATE TRANSIENT PARAMETERS Default: 1.0

  • glu_amp: Default: 1000.0

  • glu_cv: Default: 0.5

  • N_SK: SK channels Default: 15

  • SK_gamma: Default: 0.01

  • SK_Erev: Default: -90

  • SK_gating_half: Default: 0.33

  • SK_time: Default: 6.3

  • SK_hill: Default: 6

  • p_SK_bcwd: Default: (0.09391588258147192, 98.85165844770867, -147.61669527876904, 149.37767054612135)

  • bcwd_SK: Default: pSKbcwd[4] + pSKbcwd[3] / (1 + exp(pSKbcwd[1] * (temprates - pSK_bcwd[2])))

  • p_SK_frwd: Default: (-0.334167923607112, 25.590920461511878, 2.2052151559841193, 0.005904170174699533)

  • frwd_SK: Default: pSKfrwd[4] + pSKfrwd[3] / (1 + exp(pSKfrwd[1] * (temprates - pSK_frwd[2])))

  • CaM_con: **Concentrations Default: 30.0

  • mKCaM_con: Default: 70.0

  • mCaN_con: Default: 20.0

  • kon_1C: Chang Pepke model - CaM reactions I Default: 0.005

  • koff_1C: Default: 0.05

  • kon_2C: Default: 0.01

  • koff_2C: Default: 0.01

  • kon_1N: Default: 0.1

  • koff_1N: Default: 2.0

  • kon_2N: Default: 0.2

  • koff_2N: Default: 0.5

  • kf_CaM0: Chang Pepke model - CaM reactions II Default: 3.8e-6

  • kb_CaM0: Default: 0.0055

  • kf_CaM2C: Default: 0.00092

  • kb_CaM2C: Default: 0.0068

  • kf_CaM2N: Default: 0.00012

  • kb_CaM2N: Default: 0.0017

  • kf_CaM4: Default: 0.03

  • kb_CaM4: Default: 0.0015

  • kon_K1C: Chang Pepke model - CaMKII reactions Default: 0.044

  • koff_K1C: Default: 0.033

  • kon_K2C: Default: 0.044

  • koff_K2C: Default: 0.0008

  • kon_K1N: Default: 0.076

  • koff_K1N: Default: 0.3

  • kon_K2N: Default: 0.076

  • koff_K2N: Default: 0.02

  • p_camkii_q10: Chang Pepke model - autophosphorilation Default: (0.5118207068695309, 45.47503600542303, -161.42634157226917, 162.1718925882677)

  • q10: Default: pcamkiiq10[4] + pcamkiiq10[3] / (1 + exp(pcamkiiq10[1] * (temprates - pcamkii_q10[2])))

  • k1: Default: 0.0126

  • k2: Default: q10 * 0.00033

  • k3: Default: 4 * q10 * 0.00017

  • k4: Default: 4 * 4.1e-5

  • k5: Default: 4 * q10 * 2 * 1.7e-5

  • p_CaN_frwd: CaM-CaN reactions Default: (-0.29481489145354556, 29.999999999999968, 0.15940019940354327, 0.870299900298228)

  • kcanf: Default: (pCaNfrwd[4] + pCaNfrwd[3] / (1 + exp(pCaNfrwd[1] * (temprates - pCaN_frwd[2])))) * 0.0175

  • p_CaN_bcwd: Default: (-0.6833299932488973, 26.277500129849113, 0.7114524682690591, 0.29037766196937326)

  • kcanb: Default: (pCaNbcwd[4] + pCaNbcwd[3] / (1 + exp(pCaNbcwd[1] * (temprates - pCaN_bcwd[2])))) * 2.0e-5

  • p_frwd_VGCC: VGCCs Default: (1.0485098341579628, 30.66869198447378, -0.3040010721391852, 2.5032059559264357)

  • frwd_VGCC: Default: pfrwdVGCC[4] + pfrwdVGCC[3] / (1 + exp(pfrwdVGCC[1] * (temprates - pfrwd_VGCC[2])))

  • p_bcwd_VGCC: Default: (-0.3302682317933842, 36.279019647221226, 3.2259761593440155, 0.7298285671937866)

  • bcwd_VGCC: Default: pbcwdVGCC[4] + pbcwdVGCC[3] / (1 + exp(pbcwdVGCC[1] * (temprates - pbcwd_VGCC[2])))

  • Erev_CaT: Default: 10.0

  • Erev_CaR: Default: 10.0

  • Erev_CaL: Default: 10.0

  • gamma_CaT: Default: 0.012

  • gamma_CaR: Default: 0.017

  • gamma_CaL: Default: 0.027

  • N_caT: Default: 3

  • N_caR: Default: 3

  • N_caL: Default: 3

  • EGTA_kf: Calcium dye and buffers Default: 0.0027

  • EGTA_kb: Default: 0.18 * EGTA_kf

  • EGTA_con: Default: 0.0

  • BAPTA_kf: Default: 0.45

  • BAPTA_kb: Default: 0.176BAPTA_kf

  • BAPTA_con: Default: 0.0

  • Imbuf_k_on: Default: 0.247

  • Imbuf_k_off: Default: 0.524

  • K_buff_diss: Default: Imbufkoff / Imbufkon

  • Imbuf_con: Default: 62

  • Imbuf_con_dend: Default: Imbuf_con * 4

  • ogb1_kf: Calcium fluorescent dyes Default: 0.8

  • ogb1_kb: Default: 0.16

  • fluo4_kf: Default: 0.8

  • fluo4_kb: Default: 0.24

  • dye: Default: 0.0

  • fluo5f_kf: Default: dye * 0.01

  • fluo5f_kb: Default: dye * 26 * fluo5f_kf

  • fluo5f_con: Default: dye * 200.0

source

Protocols

SynapseElife.firingPatternFunction
firingPattern(
;
    start_time,
    n_pre,
    delay_pre,
    n_pos,
    delay_pos,
    delay,
    pulse,
    freq,
    causal,
    repeat_times,
    repeat_after
)

Generate external stimulation times and indices for pre/post stimulation. Usually used with dataProtocol (see folder examples in examples/).

Output

  • event_times sorted list of external event times
  • is_pre_or_post_index whether the external events are pre (true) or post (false)
source
SynapseElife.dataProtocolFunction
dataProtocol(paper)

Structure to describe a plasticity protocol "conoc"

Arguments

  • paper can be ["oconnor", "Bittner", "Goldings01", "Buchenan", "Tigaret_jitter_timespent", "Tigaret_jitter_double_1", "Poisson", "Tigaret_jitter", "Dudek_jitter", "TigaretMellor16", "TigaretMellor16_poisson", "Sleep_age", "Poisson_physiological_range", "YannisDebanne20_freq", "YannisDebanne20_freq_delay", "YannisDebanne20_ratio", "YannisDebanne20_inv", "Tigaret_burst", "Tigaret_burst_temp", "Tigaret_burst_age", "Tigaret_burst_freq", "Tigaret_burst_ca", "Tigaret_burst_dist", "Tigaret_freq_1200", "Tigaret_tripost", "Tigaret_preburst", "Tigaret_single", "TigaretMellor_sparse", "TigaretMellor_jitter_sparse", "DudekBear92-BCM-Ca", "DudekBear92-BCM", "DudekBear92-BCM-900", "DudekBear92-BCM-37", "DudekBear92-BCM-33", "DudekBear92-BCM-priming", "DudekBear92_timespend", "DudekBear92-100", "DudekBear_short", "FujiBito", "DudekBear92-sliding", "DudekBear92_temp", "DudekBear92_Ca", "DudekBear92_Mg", "DudekBear92_dist", "DudekBear92-Age", "Blocking_age_control", "Blocking_yNMDA", "Blocking_oNMDA", "Blocking_yBaP", "Blocking_oBaP", "Blocking_yGABA", "Blocking_oGABA", "DudekBear92_BCM_recovery", "DudekBear93-LFS", "DudekBear93-TBS", "Cao-TBS", "RecoverLTD", "Chang19", "Fujii_CaN", "YannisDebanne20", "YannisDebanne_temp", "YannisDebanne_age", "Meredith03-GABA", "TigaretvsMeredith", "YannisvsMeredith", "WittenbergWang06_D", "WittenbergWang06_B", "WittenbergWang06_P", "Mizuno01-LTP-Mg", "Mizuno01LTP", "Mizuno01LTD"]
  • n_pre number of presynaptic pulses
  • delay_pre delay between presynaptic pulses
  • n_pos number of postsynaptic pulses
  • delay_pos delay between postsynaptic pulses
  • delay delay between pre and postsynaptic spikes (used in STDP)
  • pulse number of pulses/pairing repetitions
  • freq frequency
  • causal causal inverts the order of pre and post, uses true or false
  • protocol name of the protocol
  • weight info entry with weight change value (not used in the model)
  • outcome String to show the qualitative outcome (not used in the model)
  • paper paper is a String to choose a predefined protocol (e.g. paper = TigaretMellor16)
  • temp temperature by which all temperature-sensitive mechanisms will be adapted
  • injection (nA) current injected in the soma, used to make postsynaptic spikes (BaPs)
  • exca (μM) extracellular calcium, we expressed it in μM to be used in the ghk function
  • exmg (mM) extracellular magnesium, we expressed it in mM, however, it is converted to μM in the magnesium relief function
  • repeat_times number of epochs additional epochs to include, for instance, TBS is usually expressed in epochs
  • repeat_after time difference between the epochs (ms)
  • testing_freq some protocols use test frequencies, this can be useful to evaluate the effect or absence of it, or to add a regular background presynaptic stimuli. 0 Hz will turn it off
  • inj_time duration (ms) of the current injection to elicit a BaP
  • age animal age (rat)
  • AP_by_EPSP yes or no, to let additional BaPS induced by EPSPs to be included in the stimulation
  • GABA_block yes or no GABAr conductance is set zero
  • jitter standard deviation of a Gaussian used to jitter the spikes
  • sparse percentage of spikes skipped
  • dista distance from the soma in μm
source

Pre-synaptic simulation

SynapseElife.stpFunction
stp(
    t_end,
    param,
    all_events_times,
    is_pre_or_post_index;
    _plot,
    nu_stp,
    kwargs...
)

This function performs the simulation of the presynaptic side.

Arguments:

  • t_end end of simulation time
  • param named tuple with parameters, example: (τ_rec = 20000, δ_decay = .0004, tau_pre = 20, tau_soma = 40)
  • all_events_times sorted list of times of external events
  • is_pre_or_post_index whether the external events are pre (true) or post (false)

Keyword arguments

  • _plot = false whether to plot the result or not
  • nu_stp = stp_build_transition_matrix() transition matrix

Output

  • is_glu_release true (success) or false (failure), list of presynaptic stimuli in all_events_index that led to a vesicle (glutamate) release
  • D, also (XD[2,:]), the evolution of the docked pool
  • R, also (XD[3,:]), the evolution of the reserve pool
  • time simulation times
  • glu_release_times the vector of times in which a successful releases (glutamate) occurred
  • bap_by_epsp_times the vector of times in which an AP was triggered by an EPSP
source

Post-synaptic simulation

SynapseElife.evolveSynapseFunction
evolveSynapse(
    xc0,
    xd0,
    p_synapse,
    events_sorted_times,
    is_pre_or_post_event,
    bap_by_epsp,
    is_glu_released,
    algos;
    verbose,
    progress,
    abstol,
    reltol,
    save_positions,
    nu,
    kwargs...
)

Perform a simulation of the synapse model. Among other things, you need to provide the external events impacting the synapse: Glutamate releases and BaPs.

Arguments

  • xc0 initial condition for the continuous variables. Example xc0 = Synapse.initial_conditions_continuous_temp(p_synapse)
  • xd0 initial condition for the discrete variables. Example xd0 = Synapse.initial_conditions_discrete(p_synapse)
  • p_synapse::SynapseParams synapse parameters. Example p_synapse = SynapseParams().
  • events_sorted_times sorted list of times (ms) for external events (Glutamate / BaP)
  • is_pre_or_post_event::Vector{Bool} whether the corresponding event in events_sorted_times is a Glutamate event. If false, it corresponds to a BaP event.
  • bap_by_epsp::Vector{<:Real} Additionnal BaPs time events, these are evoked by EPSPs. There are added to the ones indexed in is_events_glu.
  • is_glu_released::Vector{Bool} variable to shut down the Glutamate event, i.e. make the glutamate amplitude be zero. This proves useful to have this variable for reproducing some experiments. If equals to false, then glutamate amplitude is set to zero.
  • algos simulation algorithms from PiecewiseDeterministicMarkovProcesses. For example (PDMP.CHV(:lsoda), PDMP.CHV(:lsoda))

Optional arguments

  • verbose = false display information during simulation
  • abstol = 1e-8 absolute tolerance for ODE time stepper
  • reltol = 1e-7 relative tolerance for ODE time stepper
  • progress = false show a progressbar during simulation
  • save_positions = (false, true) save the values (before, after) the jumps (transitions)
  • nu transition matrix. It is initialised with buildTransitionMatrix().
source

Utils

SynapseElife.statistics_jumpsFunction

Returns the number of jumps for each component and t ∈ stpan. Usefull to determine who jumps the most. The syntax is as follows

statistics_jumps(t, xd; tspan)

where t is time, xd the discrete component and It is an inferior bound on t. You can use it like

statistics_jumps(t, out[:ampa])
source
SynapseElife.get_namesFunction

Returns a dictionary out from a simulation with associated names like out[:ampa]. Use it like

get_names(xc, xd)

Another use is for plotting.

source
SynapseElife.plot_variableFunction
plot_variable(tt, xc, xd; ...)
plot_variable(tt, xc, xd, s; tspan, kwargs...)

Plotting function

Arguments

  • tt times
  • xc continuous variable
  • xd discrete variable
  • s::Symbol = :ampa which variable to plot. Must be :ampa, :nmda, :vgcc_t, :vgcc_r, :vgcc_l or :Vsp,:Vdend,:Vsoma,:λ,:ImbufCa,:Ca,:Dye,:CaM0,:CaM2C,:CaM2N,:CaM4,:mCaN,:CaN4,:mKCaM,:KCaM0,:KCaM2N,:KCaM2C,:KCaM4,:PCaM0,:PCaM2C,:PCaM2N,:PCaM4,:P,:P2,:LTD,:LTP,:LTD_act,:LTP_act,:m,:h,:n,:SK,:λ_age,:λ_aux.

Optional arguments

  • all arguments from Plots.jl. For example xlims=(0,10), legend=false
source
SynapseElife.plot_discreteFunction
plot_discrete(tt, xc, xd; tspan, kwargs...)

Plot all discrete variables.

Arguments

  • tt times
  • xc continuous variables (result from PDMP)
  • xd discrete variables (result from PDMP)

Optional arguments

  • all arguments from Plots.jl. For example xlims = (0, 10), legend = false
source

Automatic Code generation

SynapseElife.extractReactionsCamKIIFunction
extractReactionsCamKII()
extractReactionsCamKII(filename)

This function extracts the differential equations from the CaM-CaMKII-CaN reactions and write it into the file with name filename = "write-equation.txt". The outcome of this function is used to write the F_synapse.

source