Xyce  6.1
N_DEV_Neuron6.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------------
2 // Copyright Notice
3 //
4 // Copyright 2002 Sandia Corporation. Under the terms
5 // of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S.
6 // Government retains certain rights in this software.
7 //
8 // Xyce(TM) Parallel Electrical Simulator
9 // Copyright (C) 2002-2015 Sandia Corporation
10 //
11 // This program is free software: you can redistribute it and/or modify
12 // it under the terms of the GNU General Public License as published by
13 // the Free Software Foundation, either version 3 of the License, or
14 // (at your option) any later version.
15 //
16 // This program is distributed in the hope that it will be useful,
17 // but WITHOUT ANY WARRANTY; without even the implied warranty of
18 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU General Public License
22 // along with this program. If not, see <http://www.gnu.org/licenses/>.
23 //-----------------------------------------------------------------------------
24 
25 //-------------------------------------------------------------------------
26 // Filename : $RCSfile: N_DEV_Neuron6.C,v $
27 //
28 // Purpose :
29 //
30 // Special Notes :
31 //
32 // Creator : Richard Schiek, Electrical and Microsytem Modeling
33 //
34 // Creation Date : 06/10/09
35 //
36 // Revision Information:
37 // ---------------------
38 //
39 // Revision Number: $Revision: 1.72 $
40 //
41 // Revision Date : $Date: 2015/10/20 22:32:48 $
42 //
43 // Current Owner : $Author: tvrusso $
44 //-------------------------------------------------------------------------
45 
46 #include <Xyce_config.h>
47 
48 
49 // ---------- Standard Includes ----------
50 
51 
52 // ---------- Xyce Includes ----------
53 #include <N_DEV_DeviceOptions.h>
54 #include <N_DEV_DeviceMaster.h>
55 #include <N_DEV_ExternData.h>
56 #include <N_DEV_MatrixLoadData.h>
57 #include <N_DEV_MembraneCS.h>
58 #include <N_DEV_MembraneHH.h>
59 #include <N_DEV_MembranePassive.h>
61 #include <N_DEV_Neuron6.h>
63 #include <N_DEV_SolverState.h>
64 #include <N_DEV_Message.h>
65 #include <N_ERH_ErrorMgr.h>
66 
67 #include <N_DEV_Neuron.h>
68 
69 #include <N_LAS_Vector.h>
70 #include <N_LAS_Matrix.h>
71 #include <N_UTL_FeatureTest.h>
72 #include <N_UTL_Math.h>
73 
74 namespace Xyce {
75 namespace Device {
76 
77 namespace Neuron6 {
78 
80 {
81  p.addPar ("R",1.0,&Neuron6::Instance::rInt)
82  .setGivenMember(&Neuron6::Instance::rIntGiven)
83  .setUnit(U_OHMM)
84  .setCategory(CAT_NONE)
85  .setDescription("Intracellular resistivity");
86  // typical values 1-3 kOhm-mm; use 1 Kohm-mm which is the same as 1 Ohm-m
87  p.addPar ("A",0.00025,&Neuron6::Instance::radius)
88  .setGivenMember(&Neuron6::Instance::radiusGiven)
89  .setUnit(U_METER)
90  .setCategory(CAT_NONE)
91  .setDescription("Segment radius");
92  // 250 microns,based on NEURON default diameter of 500 microns
93  p.addPar ("L",0.0001,&Neuron6::Instance::length)
94  .setGivenMember(&Neuron6::Instance::lengthGiven)
95  .setUnit(U_METER)
96  .setCategory(CAT_NONE)
97  .setDescription("Cable length");
98  // 100 microns,based on NEURON default length
99 
101  .setGivenMember(&Neuron6::Instance::nSegGiven)
102  .setUnit(U_NONE)
103  .setCategory(CAT_NONE)
104  .setDescription("Number of segments");
105 }
106 
108 {
109  p.addPar ("CMEM",0.0,&Neuron6::Model::cMem)
110  .setGivenMember(&Neuron6::Model::cMemGiven)
111  .setUnit(U_FARADMM2)
112  .setCategory(CAT_NONE)
113  .setDescription("Membrane capacitance");
114 
115  p.addPar ("GMEM",0.0,&Neuron6::Model::gMem)
116  .setGivenMember(&Neuron6::Model::gMemGiven)
117  .setUnit(U_OHMM1MM2)
118  .setCategory(CAT_NONE)
119  .setDescription("Membrane conductance");
120 
121  p.addPar ("VREST",0.0,&Neuron6::Model::vRest)
122  .setGivenMember(&Neuron6::Model::vRestGiven)
123  .setUnit(U_VOLT)
124  .setCategory(CAT_NONE)
125  .setDescription("Resting potential");
126 
127  p.addPar ("EK",0.0,&Neuron6::Model::eK)
128  .setGivenMember(&Neuron6::Model::eKGiven)
129  .setUnit(U_VOLT)
130  .setCategory(CAT_NONE)
131  .setDescription("Potassium resting potential");
132 
133  p.addPar ("GK",0.0,&Neuron6::Model::gK)
134  .setGivenMember(&Neuron6::Model::gKGiven)
135  .setUnit(U_OHMM1MM2)
136  .setCategory(CAT_NONE)
137  .setDescription("Potassium base conductance");
138 
139  p.addPar ("ENA",0.0,&Neuron6::Model::eNa)
140  .setGivenMember(&Neuron6::Model::eNaGiven)
141  .setUnit(U_VOLT)
142  .setCategory(CAT_NONE)
143  .setDescription("Sodium resting potential");
144 
145  p.addPar ("GNA",0.0,&Neuron6::Model::gNa)
146  .setGivenMember(&Neuron6::Model::gNaGiven)
147  .setUnit(U_OHMM1MM2)
148  .setCategory(CAT_NONE)
149  .setDescription("Sodium base conductance");
150 
151  p.addPar ("EA",0.0,&Neuron6::Model::eA)
152  .setGivenMember(&Neuron6::Model::eAGiven)
153  .setUnit(U_VOLT)
154  .setCategory(CAT_NONE)
155  .setDescription("a-current rest potential");
156 
157  p.addPar ("GA",0.0,&Neuron6::Model::gA)
158  .setGivenMember(&Neuron6::Model::gAGiven)
159  .setUnit(U_OHMM1MM2)
160  .setCategory(CAT_NONE)
161  .setDescription("a-current base conductance");
162 
163  p.addPar ("ECA",0.0,&Neuron6::Model::eCa)
164  .setGivenMember(&Neuron6::Model::eCaGiven)
165  .setUnit(U_VOLT)
166  .setCategory(CAT_NONE)
167  .setDescription("Calcium rest potential");
168 
169  p.addPar ("GCA",0.0,&Neuron6::Model::gCa)
170  .setGivenMember(&Neuron6::Model::gCaGiven)
171  .setUnit(U_OHMM1MM2)
172  .setCategory(CAT_NONE)
173  .setDescription("Calcium base conductance");
174 
175  p.addPar ("EKCA",0.0,&Neuron6::Model::eKCa)
176  .setGivenMember(&Neuron6::Model::eKCaGiven)
177  .setUnit(U_VOLT)
178  .setCategory(CAT_NONE)
179  .setDescription("Potassium-calcium rest potential");
180 
181  p.addPar ("GKCA",0.0,&Neuron6::Model::gKCa)
182  .setGivenMember(&Neuron6::Model::gKCaGiven)
183  .setUnit(U_OHMM1MM2)
184  .setCategory(CAT_NONE)
185  .setDescription("Potassium-calcium base conductance");
186 
187  p.addPar ("CAINIT",0.0,&Neuron6::Model::CaInit)
188  .setGivenMember(&Neuron6::Model::CaInitGiven)
189  .setUnit(U_MOLAR)
190  .setCategory(CAT_NONE)
191  .setDescription("initial intra-cellular calcium concentration");
192 
193  p.addPar ("CAGAMMA",0.0,&Neuron6::Model::CaGamma)
194  .setGivenMember(&Neuron6::Model::CaGammaGiven)
195  .setUnit(U_NONE)
196  .setCategory(CAT_NONE)
197  .setDescription("calcium current to concentration multiplier");
198 
199  p.addPar ("CATAU",0.0,&Neuron6::Model::CaTau)
200  .setGivenMember(&Neuron6::Model::CaTauGiven)
201  .setUnit(U_SECOND)
202  .setCategory(CAT_NONE)
203  .setDescription("calcium removal time constant");
204 
205  p.addPar ("R",1.0,&Neuron6::Model::rInt)
206  .setGivenMember(&Neuron6::Model::rIntGiven)
207  .setUnit(U_OHMM)
208  .setCategory(CAT_NONE)
209  .setDescription("Intracellular resistivity");
210  // typical values 1-3 kOhm-mm; use 1 Kohm-mm which is the same as 1 Ohm-m
211 
212  p.addPar ("A",0.00025,&Neuron6::Model::radius)
213  .setGivenMember(&Neuron6::Model::radiusGiven)
214  .setUnit(U_METER)
215  .setCategory(CAT_NONE)
216  .setDescription("Segment radius");
217  // 250 microns,based on NEURON default diameter of 500 microns
218 
219  p.addPar ("L",0.0001,&Neuron6::Model::length)
220  .setGivenMember(&Neuron6::Model::lengthGiven)
221  .setUnit(U_METER)
222  .setCategory(CAT_NONE)
223  .setDescription("Cable length");
224  // 100 microns,based on NEURON default length
225 
226  p.addPar ("I",0.0,&Neuron6::Model::I)
227  .setExpressionAccess(ParameterType::SOLN_DEP)
228  .setUnit(U_AMP)
229  .setCategory(CAT_NONE)
230  .setDescription("Current for user-defined current equation");
231 
232  p.addPar("IONCHANNELMODEL","",&Neuron6::Model::ionChannelModel)
233  .setGivenMember(&Neuron6::Model::ionChannelModelGiven)
234  .setUnit(U_NONE)
235  .setCategory(CAT_NONE)
236  .setDescription("Neuron6::Model to use for ion channels");
237 
238  p.addPar("N",1,&Neuron6::Model::nSeg)
239  .setGivenMember(&Neuron6::Model::nSegGiven)
240  .setUnit(U_NONE)
241  .setCategory(CAT_NONE)
242  .setDescription("Number of segments");
243 
244  p.addPar("MM_CURRENT",std::vector<std::string>(),&Neuron6::Model::membraneCurrentEqus)
246  .setUnit(U_NONE)
247  .setCategory(CAT_NONE)
248  .setDescription("Contribution to membrane current");
249 
250  p.addPar("MM_INDVARS",std::vector<std::string>(),&Neuron6::Model::membraneIndpVars)
251  .setGivenMember(&Neuron6::Model::membraneIndpVarsGiven)
252  .setUnit(U_NONE)
253  .setCategory(CAT_NONE)
254  .setDescription("Independant variables for ion channel equations");
255 
256  p.addPar("MM_INDFEQUS",std::vector<std::string>(),&Neuron6::Model::membraneIndpFEqus)
258  .setUnit(U_NONE)
259  .setCategory(CAT_NONE)
260  .setDescription("Independant variables: F equations");
261 
262  p.addPar("MM_INDQEQUS",std::vector<std::string>(),&Neuron6::Model::membraneIndpQEqus)
264  .setUnit(U_NONE)
265  .setCategory(CAT_NONE)
266  .setDescription("Independant variables: Q equations");
267 
268  p.addPar("MM_FUNCTIONS",std::vector<std::string>(),&Neuron6::Model::membraneFunctions)
270  .setUnit(U_NONE)
271  .setCategory(CAT_NONE)
272  .setDescription("Functions for membrane Neuron6::Model");
273 
274  p.addPar("MM_PARAMETERS",std::vector<std::string>(),&Neuron6::Model::membraneParameters)
276  .setUnit(U_NONE)
277  .setCategory(CAT_NONE)
278  .setDescription("Parameters for membrane Neuron6::Model");
279 }
280 
281 
282 // Class Instance
283 
284 //-----------------------------------------------------------------------------
285 // Function : Instance::Instance
286 // Purpose : constructor
287 // Special Notes :
288 // Scope : public
289 // Creator : Richard Schiek, Electrical and Microsytem Modeling
290 // Creation Date : 06/10/09
291 //-----------------------------------------------------------------------------
293  const Configuration & configuration,
294  const InstanceBlock & IB,
295  Model & Miter,
296  const FactoryBlock & factory_block)
297  : DeviceInstance(IB, configuration.getInstanceParameters(), factory_block),
298  model_(Miter),
299  rInt(0.0),
300  radius(0.0),
301  length(0.0),
302  segArea(0.0),
303  gSeg(0.0),
304  nSeg(0),
305  rIntGiven(false),
306  radiusGiven(false),
307  lengthGiven(false),
308  nSegGiven(false),
309  numIntVarsPerSegment(0),
310  numStateVarsPerSegment(0),
311  kcl1Fvalue(0.0),
312  kcl2Fvalue(0.0),
313  dkcl1F_dVin(0.0),
314  dkcl1F_dVs0(0.0),
315  dkcl2F_dVout(0.0),
316  dkcl2F_dVsn(0.0),
317  li_Pos(0),
318  li_Neg(0),
319  APosEquPosNodeOffset(0),
320  APosEquNextNodeOffset(0),
321  ANegEquNegNodeOffset(0),
322  ANegEquLastNodeOffset(0)
323 {
324  numExtVars = 2; // input and output voltage
325 
326  // Set params to constant default values:
327  setDefaultParams ();
328 
329  // Set params according to instance line
330  setParams (IB.params);
331 
332  // Set any non-constant parameter defaults:
333 
334  //if (!given("TEMP"))
335  // temp = getDeviceOptions().temp.dVal();
336 
337  // Calculate any parameters specified as expressions:
339 
340  // calculate dependent (ie computed) params and check for errors:
341  processParams ();
342 
343  // now we can calculate the number of internal vars
344  // by default we'll have a voltage at each node (nSeg)
345  // and no state vars. If the user has an ion channel on, then we'll add in vars for that
346  // ask the membrane model for the number of vars it has.
349 
350  /*
351  if( model_.ConnorStevensOn_ ) {
352  numIntVarsPerSegment += 9;
353  numStateVarsPerSegment += 2; // two currents per segment
354  }
355  */
356 
359 
360  // total up number of vars.
361  int numVars = numExtVars + numIntVars;
362 
363  //
364  // i(n) - I(n)/A - g(n,n+1) * (V(n+1) - V(n)) - g(n,n-1) * (V(n-1) - V(n)) + Cm dV(n)/d(t) = 0 LHS needs leak term memG ( vSeg - vRest )
365  //
366  // Vin : i(n) - I(n)/A - g(n,n+1) * (V(n+1) - V(n)) + Cm dV(n)/d(t) = 0
367  // Vout : i(n) - I(n)/A - g(n,n-1) * (V(n-1) - V(n)) + Cm dV(n)/d(t) = 0
368  // Vnode : i(n) - I(n)/A - g(n,n+1) * (V(n+1) - V(n)) - g(n,n-1) * (V(n-1) - V(n)) + Cm dV(n)/d(t) = 0
369  // plus node supporting equations (a, b, m)
370  //
371  // jacobian format for just the passive cable
372  // Vin Vout V1 V(nSeg)
373  // kcl Vin yes g(0,1)
374  // kcl Vout yes g(n,n-1)
375  // kcl V1 yes yes yes
376  //
377 
378  // if the membrane model includes additional internal variables, the above changes to something
379  // along these lines:
380  // Vin Vout V1 x1 y1 V(nSeg) x(nSeg) y(nSeg)
381  // kcl Vin yes g(0,1)
382  // kcl Vout yes g(n,n-1)
383  // kcl V1 yes yes yes yes yes
384  // note that V1's dependence on internal variables for segment 1 comes before dependence on
385  // Vnext, except for the last segment, in which case Vnext is Vout
386 
387  // set up jacStamp. This is dependant on the membrane model. The only part this
388  // constructor really knows about is the external variables Vin and vOut and internal node
389  // voltages
390 
391  if( jacStamp.empty() ) // redundant as jacStamp is not static for this device
392  { // it can't be as each cable may have a different number of nodes
393  jacStamp.resize(numVars);
394  jacStamp[0].resize(2);
395  jacStamp[0][0] = 0; // Vin
396  jacStamp[0][1] = 2; // Vseg[0]
397  jacStamp[1].resize(2);
398  jacStamp[1][0] = 1; // Vout
399  jacStamp[1][1] = numVars - numIntVarsPerSegment; // VnSeg[nSeg-1]
400 
401  // now loop over each segment and have the membrane model fill that instanceContainer
402  for( int i=0; i<nSeg; i++ )
403  {
404  // the cable model should take care of the Vpre, V, Vnext dependence as that
405  // is part of the cable equation. Let the membraneModel_ handle what happens
406  // at the membrane level
407  int offset = numExtVars + i * numIntVarsPerSegment; // row for vSeg
408  jacStamp[offset].resize( numIntVarsPerSegment + 2 ); // + 2 for Vin and Vout
409 
410  // have to handle a few special cases which can alter the ordering of Vprev, Vseg, Vnext
411  // for each segment number, save the offsets for the previous, current, and next segment so
412  // we don't have to rethink these special cases every time.
413  if( nSeg == 1 )
414  {
415  // here the ordering is Vin Vout Vseg
416  prevMap[i] = 0;
417  nextMap[i] = 1;
418  segMap[i] = 2;
419  jacStamp[offset][prevMap[i]] = 0; // Vin
420  jacStamp[offset][nextMap[i]] = 1; // Vout
421  jacStamp[offset][segMap[i]] = 2; // Vseg
422 
423  }
424  else if( i==0 )
425  {
426  // ordering is Vin Vseg (seg int vars) Vnext
427  prevMap[i] = 0;
428  segMap[i] = 1;
429  nextMap[i] = numIntVarsPerSegment+1;
430  jacStamp[offset][prevMap[i]] = 0; // Vin
431  jacStamp[offset][segMap[i]] = offset; // Vseg
432  jacStamp[offset][nextMap[i]] = offset + numIntVarsPerSegment; // Vnext
433  }
434  else if( i==(nSeg-1) )
435  {
436  // ordering is Vout Vprev Vseg
437  nextMap[i] = 0;
438  prevMap[i] = 1;
439  segMap[i] = 2;
440  jacStamp[offset][nextMap[i]] = 1; // Vout
441  jacStamp[offset][prevMap[i]] = offset - numIntVarsPerSegment; // Vprev
442  jacStamp[offset][segMap[i]] = offset; // Vseg
443  }
444  else
445  {
446  // ordering is Vprev Vseg (seg int vars) Vnext
447  prevMap[i] = 0;
448  segMap[i] = 1;
449  nextMap[i] = numIntVarsPerSegment+1;
450  jacStamp[offset][prevMap[i]] = offset - numIntVarsPerSegment; // Vprev
451  jacStamp[offset][segMap[i]] = offset; // Vseg
452  jacStamp[offset][nextMap[i]] = offset + numIntVarsPerSegment; // Vnext
453  }
454 
455  // pass the membraneModel_ enough information for it to construct its part of the jacobian
456  model_.membraneModel_->setJacStamp( numExtVars, i, segMap[i], jacStamp );
457  }
458 
459  }
460 
461  /*
462  // print out jacStamp
463  Xyce::dout() << "jacStamp for Neuron6" << std::endl;
464  int numRows = jacStamp.size();
465  for( int i=0; i< numRows; i++ )
466  {
467  int numCol = jacStamp[i].size();
468  Xyce::dout() << "jacStamp[ " << i << " ] = { ";
469  for(int j=0; j<numCol; j++)
470  {
471  Xyce::dout() << jacStamp[i][j] << " ";
472  }
473  Xyce::dout() << " } " << std::endl;
474  }
475  Xyce::dout() << std::endl;
476  */
477 
478 
479 
480  // calculate segment conductivities used in load calls:
481  // equation 6.30 Theoretical neuroscience: computational and mathematical modeling of neural systems, Peter Dayan and L. F. Abbot 2001
482  // g(n,n') = radius * (radius')^2 / ( rInt segLength * ( segLength * (radius')^2 + segLength' * (radius)^2 ) )
483  // this equation is in terms of current/area, in which case conductivity to the previous and next
484  // segment is not symmetric if the segment length and/or radius is not are not equal
485  // But since we work in current rather than current density, this is not a problem
486  double segLength = length / nSeg;
487  double rLong = rInt * segLength / (M_PI * radius * radius); // longitudinal resistance (ohm)
488 
489  // watch out for divide-by-0 cases; if resistivity is close to 0, just set conductance to some large value
490  // TODO: Really should be testing for anything close enough to zero to cause problems. Is there a better 'large' value to use?
491  if (rLong == 0.0)
492  {
493  gSeg = 1000.0;
494  }
495  else
496  {
497  gSeg = 1.0 / rLong;
498  }
499 
500  // variable indecies loads
502  {
503  li_Vol.resize(nSeg);
504  li_nPro.resize(nSeg);
505  li_mPro.resize(nSeg);
506  li_hPro.resize(nSeg);
507  li_aPro.resize(nSeg);
508  li_bPro.resize(nSeg);
509  li_MPro.resize(nSeg);
510  li_HPro.resize(nSeg);
511  li_cPro.resize(nSeg);
512  li_CaPro.resize(nSeg);
513  li_KCurrentState.resize(nSeg);
514  li_NaCurrentState.resize(nSeg);
515  segFvalue.resize(nSeg);
516  segQvalue.resize(nSeg);
517  segNEquFvalue.resize(nSeg);
518  segNEquQvalue.resize(nSeg);
519  segMEquFvalue.resize(nSeg);
520  segMEquQvalue.resize(nSeg);
521  segHEquFvalue.resize(nSeg);
522  segHEquQvalue.resize(nSeg);
523  segAEquFvalue.resize(nSeg);
524  segAEquQvalue.resize(nSeg);
525  segBEquFvalue.resize(nSeg);
526  segBEquQvalue.resize(nSeg);
527  segM_EquFvalue.resize(nSeg);
528  segM_EquQvalue.resize(nSeg);
529  segH_EquFvalue.resize(nSeg);
530  segH_EquQvalue.resize(nSeg);
531  segCEquFvalue.resize(nSeg);
532  segCEquQvalue.resize(nSeg);
533  segCaEquFvalue.resize(nSeg);
534  segCaEquQvalue.resize(nSeg);
535 
536  // jacobian elements
537  segF_dVp.resize(nSeg);
538  segF_dV.resize(nSeg);
539  segF_dVn.resize(nSeg);
540  segF_dn.resize(nSeg);
541  segF_dm.resize(nSeg);
542  segF_dh.resize(nSeg);
543  segF_da.resize(nSeg);
544  segF_db.resize(nSeg);
545  segF_dM.resize(nSeg);
546  segF_dH.resize(nSeg);
547  segF_dc.resize(nSeg);
548  segQ_dV.resize(nSeg);
549  dnF_dV.resize(nSeg);
550  dnF_dn.resize(nSeg);
551  dnQ_dn.resize(nSeg);
552  dmF_dV.resize(nSeg);
553  dmF_dm.resize(nSeg);
554  dmQ_dm.resize(nSeg);
555  dhF_dV.resize(nSeg);
556  dhF_dh.resize(nSeg);
557  dhQ_dh.resize(nSeg);
558  daF_dV.resize(nSeg);
559  daF_da.resize(nSeg);
560  daQ_da.resize(nSeg);
561  dbF_dV.resize(nSeg);
562  dbF_db.resize(nSeg);
563  dbQ_db.resize(nSeg);
564  dMF_dV.resize(nSeg);
565  dMF_dM.resize(nSeg);
566  dMQ_dM.resize(nSeg);
567  dHF_dV.resize(nSeg);
568  dHF_dH.resize(nSeg);
569  dHQ_dH.resize(nSeg);
570  dcF_dV.resize(nSeg);
571  dcF_dc.resize(nSeg);
572  dcF_dCa.resize(nSeg);
573  dcQ_dc.resize(nSeg);
574  dCaF_dV.resize(nSeg);
575  dCaF_dM.resize(nSeg);
576  dCaF_dH.resize(nSeg);
577  dCaF_dCa.resize(nSeg);
578  dCaQ_dCa.resize(nSeg);
579 
580  // state vars
581  potassiumCurrent.resize(nSeg);
582  sodiumCurrent.resize(nSeg);
583  }
584  else
585  {
586  /*
587  li_Vol.resize(nSeg);
588  segFvalue.resize(nSeg);
589  segQvalue.resize(nSeg);
590 
591  // jacobian elements
592  segF_dVp.resize(nSeg);
593  segF_dV.resize(nSeg);
594  segF_dVn.resize(nSeg);
595  segQ_dV.resize(nSeg);
596  */
597  }
598 
599 }
600 
601 //-----------------------------------------------------------------------------
602 // Function : Instance::~Instance
603 // Purpose : destructor
604 // Special Notes :
605 // Scope : public
606 // Creator : Richard Schiek, Electrical and Microsytem Modeling
607 // Creation Date : 06/10/09
608 //-----------------------------------------------------------------------------
610 {
611 }
612 
613 //-----------------------------------------------------------------------------
614 // Function : Instance::processParams
615 // Purpose :
616 // Special Notes :
617 // Scope : public
618 // Creator : Richard Schiek, Electrical and Microsytem Modeling
619 // Creation Date : 06/10/09
620 //-----------------------------------------------------------------------------
622 {
623  // If there are any time dependent parameters, set their values at for
624  // the current time.
625 
626  // pull unspecified params out of the model if they weren't specified here
627  if( !nSegGiven )
628  {
629  if ( model_.nSegGiven )
630  {
631  nSeg = model_.nSeg;
632  }
633  else // estimate it via lambda-d rule
634  {
635  // equations from The Neuron Book, pgs 122-123
636  // says d_lambda value of 0.1 is adequate for most purposes;
637  // use a smaller value if membrane time constant (tau_m) is shorter than 8 ms
638  // but generally uses directly specify nSeg if the default rule doesn't apply
639  double d_lambda = 0.1;
640  double f = 100.0; // frequency
641  // In equations from Neuron book, C was given in uF; we use F
642  double cMem = model_.cMem * 1.0e6;
643  // NEURON version of lambda_f equation took d in microns, other
644  // distances in cm, and returned lambda_f in microns:
645  // lambda_f = 1.0e5 * sqrt(2*radius/(4*M_PI*f*rInt*cMem));
646  // nSeg = int((length/(d_lambda*lambda_f)+0.9)/2)*2 + 1;
647  // I modified the coefficient in the lambda_f equation
648  // to keep distance units in cm
649  // (should also work in other units as long as the units are consistent)
650  double lambda_f = 1.0e3 * sqrt(2*radius/(4*M_PI*f*rInt*cMem));
651  nSeg = int((length/(d_lambda*lambda_f)+0.9)/2)*2 + 1;
652  }
653  }
654 
655  if( !rIntGiven && model_.rIntGiven )
656  {
657  rInt = model_.rInt;
658  }
659  if( !radiusGiven && model_.radiusGiven )
660  {
661  radius = model_.radius;
662  }
663  if( !lengthGiven && model_.lengthGiven )
664  {
665  length = model_.length;
666  }
667 
668  // now that nSeg, length and radius are set calculate segment area
669  segArea = 2.0 * M_PI * radius * length / nSeg;
670  return true;
671 }
672 
673 //-----------------------------------------------------------------------------
674 // Function : Instance::updateTemperature
675 // Purpose :
676 // Special Notes :
677 // Scope : public
678 // Creator : Richard Schiek, Electrical and Microsytem Modeling
679 // Creation Date : 06/10/09
680 //-----------------------------------------------------------------------------
681 bool Instance::updateTemperature ( const double & temp)
682 {
683  bool bsuccess = true;
684  return bsuccess;
685 }
686 
687 //-----------------------------------------------------------------------------
688 // Function : Instance::registerLIDs
689 // Purpose :
690 // Special Notes :
691 // Scope : public
692 // Creator : Richard Schiek, Electrical and Microsytem Modeling
693 // Creation Date : 06/10/09
694 //-----------------------------------------------------------------------------
695 void Instance::registerLIDs(const std::vector<int> & intLIDVecRef,
696  const std::vector<int> & extLIDVecRef)
697 {
698  AssertLIDs(intLIDVecRef.size() == numIntVars);
699  AssertLIDs(extLIDVecRef.size() == numExtVars);
700 
701  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
702  {
703  Xyce::dout() << std::endl << section_divider << std::endl;
704  Xyce::dout() << " Instance::registerLIDs" << std::endl;
705  Xyce::dout() << " name = " << getName() << std::endl;
706  }
707 
708  // copy over the global ID lists.
709  intLIDVec = intLIDVecRef;
710  extLIDVec = extLIDVecRef;
711 
712  li_Pos = extLIDVec[0];
713  li_Neg = extLIDVec[1];
714 
715  // resize our storage location for the internal vars.
716  li_internalVars.resize( numIntVars );
717 
718  // now copy in the local ids
719  for( int i=0; i<numIntVars; i++ )
720  {
721  li_internalVars[i] = intLIDVec[i];
722  }
723  /*
724  if( model_.ConnorStevensOn_ )
725  {
726  for( int i=0, j=0; i<nSeg; i++, j+=10)
727  {
728  li_Vol[i] = intLIDVec[j];
729  li_nPro[i] = intLIDVec[j+1];
730  li_mPro[i] = intLIDVec[j+2];
731  li_hPro[i] = intLIDVec[j+3];
732  li_aPro[i] = intLIDVec[j+4];
733  li_bPro[i] = intLIDVec[j+5];
734  li_MPro[i] = intLIDVec[j+6];
735  li_HPro[i] = intLIDVec[j+7];
736  li_cPro[i] = intLIDVec[j+8];
737  li_CaPro[i] = intLIDVec[j+9];
738  }
739  }
740  else
741  {
742  for( int i=0, j=0; i<nSeg; i++, j+=1)
743  {
744  li_Vol[i] = intLIDVec[j];
745  }
746  }
747 
748  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) )
749  {
750  Xyce::dout() << " li_Pos = " << li_Pos << std::endl
751  << " li_Neg = " << li_Neg << std::endl;
752  for( int i=0; i<nSeg; i++ )
753  {
754  Xyce::dout() << " li_Vol[ " << i << " ] = " << li_Vol[i] << std::endl;
755 
756  if( model_.ConnorStevensOn_ )
757  {
758  Xyce::dout() << " li_nPro[ " << i << " ] = " << li_nPro[i] << std::endl
759  << " li_mPro[ " << i << " ] = " << li_mPro[i] << std::endl
760  << " li_hPro[ " << i << " ] = " << li_hPro[i] << std::endl
761  << " li_aPro[ " << i << " ] = " << li_aPro[i] << std::endl
762  << " li_bPro[ " << i << " ] = " << li_bPro[i] << std::endl
763  << " li_MPro[ " << i << " ] = " << li_MPro[i] << std::endl
764  << " li_HPro[ " << i << " ] = " << li_HPro[i] << std::endl
765  << " li_cPro[ " << i << " ] = " << li_cPro[i] << std::endl
766  << " li_CaPro[ " << i << " ] = " << li_CaPro[i] << std::endl;
767  }
768  }
769  }
770  #endif
771 
772  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) )
773  {
774  Xyce::dout() << section_divider << std::endl;
775  }
776  #endif
777  */
778 }
779 
780 //-----------------------------------------------------------------------------
781 // Function : Instance::loadNodeSymbols
782 // Purpose :
783 // Special Notes :
784 // Scope : public
785 // Creator : Eric R. Keiter, SNL, Parallel Computational Sciences
786 // Creation Date : 05/13/05
787 //-----------------------------------------------------------------------------
788 void Instance::loadNodeSymbols(Util::SymbolTable &symbol_table) const
789 {
790  for (int i = 0; i < nSeg ; i++)
791  {
792  std::ostringstream segNumber;
793  segNumber << i;
794  addInternalNode(symbol_table, li_internalVars[i* numIntVarsPerSegment ], getName(), "V" + segNumber.str());
795  if( numIntVarsPerSegment > 1 )
796  {
797  // need to handle this better as the HH and CS (and userdefined) optinos will have differnet internal variables
798  addInternalNode(symbol_table, li_internalVars[i*numIntVarsPerSegment +1], getName(), "N" + segNumber.str());
799  addInternalNode(symbol_table, li_internalVars[i*numIntVarsPerSegment +2], getName(), "M" + segNumber.str());
800  addInternalNode(symbol_table, li_internalVars[i*numIntVarsPerSegment +3], getName(), "H" + segNumber.str());
801  }
802 
804  {
805  addInternalNode(symbol_table, li_nPro[i], getName(), "N" + segNumber.str());
806  addInternalNode(symbol_table, li_mPro[i], getName(), "M" + segNumber.str());
807  addInternalNode(symbol_table, li_hPro[i], getName(), "H" + segNumber.str());
808  addInternalNode(symbol_table, li_aPro[i], getName(), "A" + segNumber.str());
809  addInternalNode(symbol_table, li_bPro[i], getName(), "B" + segNumber.str());
810  addInternalNode(symbol_table, li_MPro[i], getName(), "M_" + segNumber.str());
811  addInternalNode(symbol_table, li_HPro[i], getName(), "H_" + segNumber.str());
812  addInternalNode(symbol_table, li_cPro[i], getName(), "C" + segNumber.str());
813  addInternalNode(symbol_table, li_CaPro[i], getName(), "Ca" + segNumber.str());
814  }
815  }
816 }
817 
818 //-----------------------------------------------------------------------------
819 // Function : Instance::registerStateLIDs
820 // Purpose :
821 // Special Notes :
822 // Scope : public
823 // Creator : Richard Schiek, Electrical and Microsytem Modeling
824 // Creation Date : 06/10/09
825 //-----------------------------------------------------------------------------
826 void Instance::registerStateLIDs( const std::vector<int> & staLIDVecRef )
827 {
828  AssertLIDs(staLIDVecRef.size() == numStateVars);
829 
830  // copy over the global ID lists.
831  staLIDVec = staLIDVecRef;
832 
834  {
835  for( int i=0, j=0; i<nSeg; i++, j+=2)
836  {
837  li_KCurrentState[i] = staLIDVec[j];
838  li_NaCurrentState[i] = staLIDVec[j+1];
839  }
840  }
841  else
842  {
843  }
844 }
845 
846 //-----------------------------------------------------------------------------
847 // Function : Instance::jacobianStamp
848 // Purpose :
849 // Special Notes :
850 // Scope : public
851 // Creator : Richard Schiek, Electrical and Microsytem Modeling
852 // Creation Date : 06/10/09
853 //-----------------------------------------------------------------------------
854 const std::vector< std::vector<int> > & Instance::jacobianStamp() const
855 {
856  return jacStamp;
857 }
858 
859 //-----------------------------------------------------------------------------
860 // Function : Instance::registerJacLIDs
861 // Purpose :
862 // Special Notes :
863 // Scope : public
864 // Creator : Richard Schiek, Electrical and Microsytem Modeling
865 // Creation Date : 06/10/09
866 //-----------------------------------------------------------------------------
867 void Instance::registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec )
868 {
869  DeviceInstance::registerJacLIDs( jacLIDVec );
870 
871  // resize our storage location and store the results
872 
873  int numRows = jacLIDVec.size();
874  jacobianOffsets.resize( numRows );
875  for( int i=0; i< numRows; i++ )
876  {
877  int numCol = jacLIDVec[i].size();
878  jacobianOffsets[i].resize( numCol );
879  for( int j=0; j< numCol; j++ )
880  {
881  jacobianOffsets[i][j] = jacLIDVec[i][j];
882  }
883  }
884 
885  // external terminals
886  APosEquPosNodeOffset = jacLIDVec[0][0];
887  APosEquNextNodeOffset = jacLIDVec[0][1];
888  ANegEquNegNodeOffset = jacLIDVec[1][0];
889  ANegEquLastNodeOffset = jacLIDVec[1][1];
890  /*
891  Xyce::dout() << "APosEquPosNodeOffset = " << APosEquPosNodeOffset << std::endl;
892  Xyce::dout() << "APosEquNextNodeOffset = " << APosEquNextNodeOffset << std::endl;
893  Xyce::dout() << "ANegEquNegNodeOffset = " << ANegEquNegNodeOffset << std::endl;
894  Xyce::dout() << "ANegEquLastNodeOffset = " << ANegEquLastNodeOffset << std::endl;
895  */
896 
898  {
899  // internal variables
900  SegVEqnVpreOffset.resize(nSeg);
901  SegVEqnVsegOffset.resize(nSeg);
902  SegVEqnVnexOffset.resize(nSeg);
903  SegVEqnNOffset.resize(nSeg);
904  SegVEqnMOffset.resize(nSeg);
905  SegVEqnHOffset.resize(nSeg);
906  NEquVNodeOffset.resize(nSeg);
907  NEquNNodeOffset.resize(nSeg);
908  MEquVNodeOffset.resize(nSeg);
909  MEquMNodeOffset.resize(nSeg);
910  HEquVNodeOffset.resize(nSeg);
911  HEquHNodeOffset.resize(nSeg);
912  AEquVNodeOffset.resize(nSeg);
913  AEquANodeOffset.resize(nSeg);
914  BEquVNodeOffset.resize(nSeg);
915  BEquBNodeOffset.resize(nSeg);
916  M_EquVNodeOffset.resize(nSeg);
917  M_EquM_NodeOffset.resize(nSeg);
918  H_EquVNodeOffset.resize(nSeg);
919  H_EquH_NodeOffset.resize(nSeg);
920  CEquVNodeOffset.resize(nSeg);
921  CEquCNodeOffset.resize(nSeg);
922  CEquCaNodeOffset.resize(nSeg);
923  CaEquVNodeOffset.resize(nSeg);
924  CaEquM_NodeOffset.resize(nSeg);
925  CaEquH_NodeOffset.resize(nSeg);
926  CaEquCaNodeOffset.resize(nSeg);
927 
928  for(int i=0, j=2; i<nSeg; i++, j+=10 )
929  {
930  SegVEqnVpreOffset[i] = jacLIDVec[j][0];
931  SegVEqnVsegOffset[i] = jacLIDVec[j][1];
932  SegVEqnNOffset[i] = jacLIDVec[j][2];
933  SegVEqnMOffset[i] = jacLIDVec[j][3];
934  SegVEqnHOffset[i] = jacLIDVec[j][4];
935  SegVEqnVnexOffset[i] = jacLIDVec[j][5];
936 
937  NEquVNodeOffset[i] = jacLIDVec[j+1][0];
938  NEquNNodeOffset[i] = jacLIDVec[j+1][1];
939  MEquVNodeOffset[i] = jacLIDVec[j+2][0];
940  MEquMNodeOffset[i] = jacLIDVec[j+2][1];
941  HEquVNodeOffset[i] = jacLIDVec[j+3][0];
942  HEquHNodeOffset[i] = jacLIDVec[j+3][1];
943  AEquVNodeOffset[i] = jacLIDVec[j+4][0];
944  AEquANodeOffset[i] = jacLIDVec[j+4][1];
945  BEquVNodeOffset[i] = jacLIDVec[j+5][0];
946  BEquBNodeOffset[i] = jacLIDVec[j+5][1];
947  M_EquVNodeOffset[i] = jacLIDVec[j+6][0];
948  M_EquM_NodeOffset[i] = jacLIDVec[j+6][1];
949  H_EquVNodeOffset[i] = jacLIDVec[j+7][0];
950  H_EquH_NodeOffset[i] = jacLIDVec[j+7][1];
951  CEquVNodeOffset[i] = jacLIDVec[j+8][0];
952  CEquCNodeOffset[i] = jacLIDVec[j+8][1];
953  CEquCaNodeOffset[i] = jacLIDVec[j+8][2];
954  CaEquVNodeOffset[i] = jacLIDVec[j+9][0];
955  CaEquM_NodeOffset[i] = jacLIDVec[j+9][1];
956  CaEquH_NodeOffset[i] = jacLIDVec[j+9][2];
957  CaEquCaNodeOffset[i] = jacLIDVec[j+9][3];
958 
959  /*
960  Xyce::dout() << SegVEqnVpreOffset[i] << ", "
961  << SegVEqnVsegOffset[i] << ", "
962  << SegVEqnNOffset[i] << ", "
963  << SegVEqnMOffset[i] << ", "
964  << SegVEqnHOffset[i] << ", "
965  << SegVEqnVnexOffset[i] << ", "
966  << NEquVNodeOffset[i] << ", "
967  << NEquNNodeOffset[i] << ", "
968  << MEquVNodeOffset[i] << ", "
969  << MEquMNodeOffset[i] << ", "
970  << HEquVNodeOffset[i] << ", "
971  << HEquHNodeOffset[i] << std::endl;
972  */
973  }
974  }
975  else
976  {
977  // internal variables
978  SegVEqnVpreOffset.resize(nSeg);
979  SegVEqnVsegOffset.resize(nSeg);
980  SegVEqnVnexOffset.resize(nSeg);
981 
982  for(int i=0, j=2; i<nSeg; i++, j+=numIntVarsPerSegment )
983  {
984  // Xyce::dout() << " i = " << i << " j = " << j << " jacLIDVec[ " << j << " ].size() = " << jacLIDVec[j].size() << std::endl;
985  SegVEqnVpreOffset[i] = jacLIDVec[j][0];
986  SegVEqnVsegOffset[i] = jacLIDVec[j][1];
987  SegVEqnVnexOffset[i] = jacLIDVec[j][2];
988 
989  /*
990  Xyce::dout() << SegVEqnVpreOffset[i] << ", "
991  << SegVEqnVsegOffset[i] << ", "
992  << SegVEqnVnexOffset[i] << std::endl;
993  */
994  }
995  }
996 
997 }
998 
999 
1000 //-----------------------------------------------------------------------------
1001 // Function : Instance::updateIntermediateVars
1002 // Purpose :
1003 // Special Notes :
1004 // Scope : public
1005 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1006 // Creation Date : 06/10/09
1007 //-----------------------------------------------------------------------------
1009 {
1010  bool bsuccess = true;
1011 
1012  return bsuccess;
1013 }
1014 //-----------------------------------------------------------------------------
1015 // Function : Instance::updatePrimaryState
1016 // Purpose :
1017 // Special Notes :
1018 // Scope : public
1019 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1020 // Creation Date : 06/10/09
1021 //-----------------------------------------------------------------------------
1023 {
1024  bool bsuccess = true;
1025 
1026  return bsuccess;
1027 }
1028 
1029 //-----------------------------------------------------------------------------
1030 // Function : Instance::updateSecondaryState
1031 // Purpose :
1032 // Special Notes :
1033 // Scope : public
1034 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1035 // Creation Date : 06/10/09
1036 //-----------------------------------------------------------------------------
1038 {
1039  bool bsuccess = true;
1040 
1041  return bsuccess;
1042 }
1043 
1044 //-----------------------------------------------------------------------------
1045 // Function : Instance::loadDAEQVector
1046 //
1047 // Purpose : Loads the Q-vector contributions for a single
1048 // Neuron 6 instance.
1049 //
1050 // Special Notes : The "Q" vector is part of a standard DAE formalism in
1051 // which the system of equations is represented as:
1052 //
1053 // f(x) = dQ(x)/dt + F(x) - B(t) = 0
1054 //
1055 // Scope : public
1056 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1057 // Creation Date : 06/10/09
1058 //-----------------------------------------------------------------------------
1060 {
1061  bool bsuccess = true;
1062 
1063  Linear::Vector * solVectorPtr = extData.nextSolVectorPtr;
1064  Linear::Vector * daeQVecPtr = extData.daeQVectorPtr;
1065 
1066  // no Q component for the cable component of this devcie
1067 
1068  // now let the membrane model load it's Q component (memCap dV/dt, etc.)
1069  for( int i=0; i<nSeg ; i++)
1070  {
1071  model_.membraneModel_->loadDAEQVector ( i, li_internalVars, solVectorPtr, daeQVecPtr, segArea );
1072  /*
1073  (*daeQVecPtr)[li_Vol[i]] += segQvalue[i];
1074  if( model_.ConnorStevensOn_ )
1075  {
1076  (*daeQVecPtr)[li_nPro[i]] += segNEquQvalue[i];
1077  (*daeQVecPtr)[li_mPro[i]] += segMEquQvalue[i];
1078  (*daeQVecPtr)[li_hPro[i]] += segHEquQvalue[i];
1079  (*daeQVecPtr)[li_aPro[i]] += segAEquQvalue[i];
1080  (*daeQVecPtr)[li_bPro[i]] += segBEquQvalue[i];
1081  (*daeQVecPtr)[li_MPro[i]] += segM_EquQvalue[i];
1082  (*daeQVecPtr)[li_HPro[i]] += segH_EquQvalue[i];
1083  (*daeQVecPtr)[li_cPro[i]] += segCEquQvalue[i];
1084  (*daeQVecPtr)[li_CaPro[i]] += segCaEquQvalue[i];
1085  }
1086  */
1087  }
1088 
1089  return bsuccess;
1090 }
1091 
1092 
1093 
1094 //-----------------------------------------------------------------------------
1095 // Function : Instance::loadDAEFVector
1096 //
1097 // Purpose : Loads the F-vector contributions for a single
1098 // Neuron 6 instance.
1099 //
1100 // Special Notes :
1101 //
1102 // Scope : public
1103 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1104 // Creation Date : 06/10/09
1105 //-----------------------------------------------------------------------------
1107 {
1108  bool bsuccess=true;
1109 
1110  Linear::Vector * solVectorPtr = extData.nextSolVectorPtr;
1111  Linear::Vector * daeFVecPtr = extData.daeFVectorPtr;
1112 
1113  double vIn = (*solVectorPtr)[li_Pos];
1114  double vOut = (*solVectorPtr)[li_Neg];
1115 
1116  // take care of the input and output nodes as they are different
1117  (*daeFVecPtr)[li_Pos] += -2.0 * gSeg * ((*solVectorPtr)[li_internalVars[0]] - vIn );
1118  (*daeFVecPtr)[li_Neg] += -2.0 * gSeg * ((*solVectorPtr)[li_internalVars[(nSeg-1)*numIntVarsPerSegment]] - vOut );
1119 
1120  for( int i=0; i<nSeg ; i++)
1121  {
1122  // for this segment get the values of the local vars
1123  double vSeg = (*solVectorPtr)[li_internalVars[i*numIntVarsPerSegment]];
1124  double vNext = 0.0;
1125  double gNext = 0.0;
1126  if (i == (nSeg - 1))
1127  {
1128  vNext = vOut;
1129  gNext = gSeg * 2.0;
1130  }
1131  else
1132  {
1133  vNext = (*solVectorPtr)[li_internalVars[(i+1)*numIntVarsPerSegment]];
1134  gNext = gSeg;
1135  }
1136  double vPrev = 0.0;
1137  double gPrev = 0.0;
1138  if (i == 0 )
1139  {
1140  vPrev = vIn;
1141  gPrev = gSeg * 2.0;
1142  }
1143  else
1144  {
1145  vPrev = (*solVectorPtr)[li_internalVars[(i-1)*numIntVarsPerSegment]];
1146  gPrev = gSeg;
1147  }
1148 
1149  // li_internalVars lists numIntVarsPerSegment for each segment. V for each segment will
1150  // be first. So, V's offset in li_internalVars is i * numIntVarsPerSegment
1151 
1152  (*daeFVecPtr)[li_internalVars[i*numIntVarsPerSegment]] += - gPrev * (vPrev - vSeg) - gNext * (vNext - vSeg);
1153 
1154  model_.membraneModel_->loadDAEFVector ( i, li_internalVars, solVectorPtr, daeFVecPtr, segArea );
1155 
1156  }
1157 
1158  return bsuccess;
1159 }
1160 
1161 //-----------------------------------------------------------------------------
1162 // Function : Instance::loadDAEdQdx
1163 //
1164 // Purpose : Loads the Q-vector contributions for a single
1165 // Neuron 6 instance.
1166 // Scope : public
1167 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1168 // Creation Date : 06/10/09
1169 //-----------------------------------------------------------------------------
1171 {
1172  bool bsuccess = true;
1173 
1174  Linear::Vector * solVectorPtr = extData.nextSolVectorPtr;
1175  Linear::Matrix * dQdxMatPtr = extData.dQdxMatrixPtr;
1176 
1177  for( int i=0; i<nSeg ; i++)
1178  {
1179  // let the membrane model load it's part
1180  model_.membraneModel_->loadDAEdQdx ( i, segMap[i], li_internalVars, jacobianOffsets, solVectorPtr, dQdxMatPtr, segArea );
1181 
1182  /*
1183  if( model_.ConnorStevensOn_ )
1184  {
1185  (*dQdxMatPtr)[li_nPro[i]][NEquNNodeOffset[i]] += dnQ_dn[i];
1186  (*dQdxMatPtr)[li_mPro[i]][MEquMNodeOffset[i]] += dmQ_dm[i];
1187  (*dQdxMatPtr)[li_hPro[i]][HEquHNodeOffset[i]] += dhQ_dh[i];
1188  (*dQdxMatPtr)[li_aPro[i]][AEquANodeOffset[i]] += daQ_da[i];
1189  (*dQdxMatPtr)[li_bPro[i]][BEquBNodeOffset[i]] += dbQ_db[i];
1190  (*dQdxMatPtr)[li_MPro[i]][M_EquM_NodeOffset[i]] += dMQ_dM[i];
1191  (*dQdxMatPtr)[li_HPro[i]][H_EquH_NodeOffset[i]] += dHQ_dH[i];
1192  (*dQdxMatPtr)[li_cPro[i]][CEquCNodeOffset[i]] += dcQ_dc[i];
1193  (*dQdxMatPtr)[li_CaPro[i]][CaEquCaNodeOffset[i]] += dCaQ_dCa[i];
1194  }
1195  */
1196  }
1197 
1198  return bsuccess;
1199 }
1200 
1201 
1202 
1203 //-----------------------------------------------------------------------------
1204 // Function : Instance::loadDAEdFdx ()
1205 //
1206 // Purpose : Loads the F-vector contributions for a single
1207 // Neuron 6 instance.
1208 //
1209 // Special Notes : This is an algebraic constaint.
1210 //
1211 // Scope : public
1212 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1213 // Creation Date : 06/10/09
1214 //-----------------------------------------------------------------------------
1216 {
1217  bool bsuccess = true;
1218 
1219  Linear::Vector * solVectorPtr = extData.nextSolVectorPtr;
1220  Linear::Matrix * dFdxMatPtr = extData.dFdxMatrixPtr;
1221 
1222  (*dFdxMatPtr)[li_Pos][APosEquPosNodeOffset] += 2.0 * gSeg;
1223  (*dFdxMatPtr)[li_Pos][APosEquNextNodeOffset] += -2.0 * gSeg;
1224 
1225  (*dFdxMatPtr)[li_Neg][ANegEquNegNodeOffset] += 2.0 * gSeg;
1226  (*dFdxMatPtr)[li_Neg][ANegEquLastNodeOffset] += -2.0 * gSeg;
1227 
1228  for( int i=0; i<nSeg ; i++)
1229  {
1230  int offset = i * numIntVarsPerSegment;
1231 
1232  int row = numExtVars + i * numIntVarsPerSegment;
1233 
1234  double gPrev = gSeg;
1235  double gNext = gSeg;
1236  if (i == 0)
1237  {
1238  gPrev = gSeg * 2.0;
1239  }
1240  if (i == (nSeg-1))
1241  {
1242  gNext = gSeg * 2.0;
1243  }
1244 
1245  (*dFdxMatPtr)[li_internalVars[offset]][jacobianOffsets[row][prevMap[i]]] += -gPrev; // Vpre
1246  (*dFdxMatPtr)[li_internalVars[offset]][jacobianOffsets[row][segMap[i]]] += gPrev + gNext; // Vseg
1247  (*dFdxMatPtr)[li_internalVars[offset]][jacobianOffsets[row][nextMap[i]]] += -gNext; // Vnext
1248 
1249  // now let the membrane model load it's part
1250  model_.membraneModel_->loadDAEdFdx ( i, segMap[i], li_internalVars, jacobianOffsets, solVectorPtr, dFdxMatPtr, segArea );
1251 
1252  }
1253 
1254  return bsuccess;
1255 }
1256 
1257 //-----------------------------------------------------------------------------
1258 // Function : Instance::setIC
1259 // Purpose :
1260 // Special Notes :
1261 // Scope : public
1262 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1263 // Creation Date : 06/10/09
1264 //-----------------------------------------------------------------------------
1266 {
1267  bool bsuccess = true;
1268 
1269  return bsuccess;
1270 }
1271 
1272 //-----------------------------------------------------------------------------
1273 // Function : Instance::varTypes
1274 // Purpose :
1275 // Special Notes :
1276 // Scope : public
1277 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1278 // Creation Date : 06/10/09
1279 //-----------------------------------------------------------------------------
1280 void Instance::varTypes( std::vector<char> & varTypeVec )
1281 {
1282  //varTypeVec.resize(1);
1283  //varTypeVec[0] = 'I';
1284 }
1285 
1286 
1287 //-----------------------------------------------------------------------------
1288 // Function : Model::Model
1289 // Purpose : block constructor
1290 // Special Notes :
1291 // Scope : public
1292 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1293 // Creation Date : 06/10/09
1294 //-----------------------------------------------------------------------------
1296  const Configuration & configuration,
1297  const ModelBlock & MB,
1298  const FactoryBlock & factory_block)
1299  : DeviceModel(MB, configuration.getModelParameters(), factory_block),
1300  cMem(0.0),
1301  gMem(0.0),
1302  vRest(0.0),
1303  eNa(0.0),
1304  gNa(0.0),
1305  eK(0.0),
1306  gK(0.0),
1307  eA(0.0),
1308  gA(0.0),
1309  eCa(0.0),
1310  gCa(0.0),
1311  eKCa(0.0),
1312  gKCa(0.0),
1313  CaInit(0.0),
1314  CaGamma(0.0),
1315  CaTau(0.0),
1316  rInt(0.0),
1317  radius(0.0),
1318  length(0.0),
1319  nSeg(0.0),
1320  rIntGiven(false),
1321  radiusGiven(false),
1322  lengthGiven(false),
1323  nSegGiven(false),
1324  ionChannelModelGiven(false),
1325  cMemGiven(false),
1326  gMemGiven(false),
1327  vRestGiven(false),
1328  eNaGiven(false),
1329  gNaGiven(false),
1330  eKGiven(false),
1331  gKGiven(false),
1332  eAGiven(false),
1333  gAGiven(false),
1334  eCaGiven(false),
1335  gCaGiven(false),
1336  eKCaGiven(false),
1337  gKCaGiven(false),
1338  CaInitGiven(false),
1339  CaGammaGiven(false),
1340  CaTauGiven(false),
1341  hodgenHuxleyOn_(false),
1342  ConnorStevensOn_(false),
1343  sodiumOn_(false),
1344  potassiumOn_(false),
1345  aCurrentOn_(false),
1346  calciumOn_(false),
1347  membraneIndpVarsGiven(false),
1348  membraneIndpFEqusGiven(false),
1349  membraneIndpQEqusGiven(false)
1350 {
1351 
1352  // Set params to constant default values:
1353  setDefaultParams ();
1354 
1355  // Set params according to .model line and constant defaults from metadata:
1356  setModParams (MB.params);
1357 
1358  // Set any non-constant parameter defaults:
1359  //if (!given("TNOM"))
1360  // tnom = getDeviceOptions().tnom;
1361 
1362  // Calculate any parameters specified as expressions:
1364 
1365  // calculate dependent (ie computed) params and check for errors:
1366 
1367  processParams ();
1368 
1369  // check the specified model type and allocate the needed membrane model
1370  if( ionChannelModelGiven )
1371  {
1372  // should change the case on ionChannelModel and simplify this set of clauses
1373  if( ionChannelModel == "passive" || ionChannelModel == "PASSIVE" )
1374  {
1376 
1377  }
1378  else if( ionChannelModel == "hh" || ionChannelModel == "HH" )
1379  {
1380  membraneModel_ = rcp( new MembraneHH( getSolverState(), cMem, gMem, vRest, eK, gK, eNa, gNa) );
1381  }
1382  else if( ionChannelModel == "cs" || ionChannelModel == "CS")
1383  {
1384  membraneModel_ = rcp( new MembraneCS( getSolverState() ) );
1385  }
1386  else if( ionChannelModel == "ud" || ionChannelModel == "UD")
1387  {
1390  }
1391  else
1392  {
1393  UserFatal0(*this) << "Model: unknown ion channel model given \"" << ionChannelModel << "\"";
1394  }
1395 
1396  }
1397  else
1398  {
1399  // no model given so assume passive
1401  }
1402 }
1403 
1404 
1405 //-----------------------------------------------------------------------------
1406 // Function : Model::~Model
1407 // Purpose : destructor
1408 // Special Notes :
1409 // Scope : public
1410 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1411 // Creation Date : 06/10/09
1412 //-----------------------------------------------------------------------------
1414 {
1415  std::vector<Instance*>::iterator iter;
1416  std::vector<Instance*>::iterator first = instanceContainer.begin();
1417  std::vector<Instance*>::iterator last = instanceContainer.end();
1418 
1419  for (iter=first; iter!=last; ++iter)
1420  {
1421  delete (*iter);
1422  }
1423 
1424 }
1425 
1426 // additional Declarations
1427 //-----------------------------------------------------------------------------
1428 // Function : Model::processParams
1429 // Purpose :
1430 // Special Notes :
1431 // Scope : public
1432 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1433 // Creation Date : 06/10/09
1434 //-----------------------------------------------------------------------------
1436 {
1437  return true;
1438 }
1439 
1440 //----------------------------------------------------------------------------
1441 // Function : Model::processInstanceParams
1442 // Purpose :
1443 // Special Notes :
1444 // Scope : public
1445 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1446 // Creation Date : 06/10/09
1447 //----------------------------------------------------------------------------
1449 {
1450 
1451  std::vector<Instance*>::iterator iter;
1452  std::vector<Instance*>::iterator first = instanceContainer.begin();
1453  std::vector<Instance*>::iterator last = instanceContainer.end();
1454 
1455  for (iter=first; iter!=last; ++iter)
1456  {
1457  (*iter)->processParams();
1458  }
1459 
1460  return true;
1461 }
1462 //-----------------------------------------------------------------------------
1463 // Function : Model::printOutInstances
1464 // Purpose : debugging tool.
1465 // Special Notes :
1466 // Scope : public
1467 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1468 // Creation Date : 06/10/09
1469 //-----------------------------------------------------------------------------
1470 std::ostream &Model::printOutInstances(std::ostream &os) const
1471 {
1472  std::vector<Instance*>::const_iterator iter;
1473  std::vector<Instance*>::const_iterator first = instanceContainer.begin();
1474  std::vector<Instance*>::const_iterator last = instanceContainer.end();
1475 
1476  int i, isize;
1477  isize = instanceContainer.size();
1478 
1479  os << std::endl;
1480  os << "Number of Neuron instances: " << isize << std::endl;
1481  os << " name=\t\tmodelName\tParameters" << std::endl;
1482  for (i=0, iter=first; iter!=last; ++iter, ++i)
1483  {
1484  os << " " << i << ": " << (*iter)->getName() << "\t";
1485  os << getName();
1486  os << std::endl;
1487  }
1488 
1489  os << std::endl;
1490  return os;
1491 }
1492 
1493 //-----------------------------------------------------------------------------
1494 // Function : Model::forEachInstance
1495 // Purpose :
1496 // Special Notes :
1497 // Scope : public
1498 // Creator : David Baur
1499 // Creation Date : 2/4/2014
1500 //-----------------------------------------------------------------------------
1501 /// Apply a device instance "op" to all instances associated with this
1502 /// model
1503 ///
1504 /// @param[in] op Operator to apply to all instances.
1505 ///
1506 ///
1507 void Model::forEachInstance(DeviceInstanceOp &op) const /* override */
1508 {
1509  for (std::vector<Instance *>::const_iterator it = instanceContainer.begin(); it != instanceContainer.end(); ++it)
1510  op(*it);
1511 }
1512 
1513 
1514 Device *Traits::factory(const Configuration &configuration, const FactoryBlock &factory_block)
1515 {
1516 
1517  return new DeviceMaster<Traits>(configuration, factory_block, factory_block.solverState_, factory_block.deviceOptions_);
1518 }
1519 
1521 {
1523  .registerDevice("neuron", 6)
1524  .registerModelType("neuron", 6);
1525 }
1526 
1527 } // namespace Neuron6
1528 } // namespace Device
1529 } // namespace Xyce
const InstanceName & getName() const
std::vector< double > segF_dM
std::map< int, int > nextMap
virtual std::ostream & printOutInstances(std::ostream &os) const
std::vector< int > SegVEqnNOffset
std::vector< double > segQ_dV
std::vector< double > dCaF_dH
std::vector< double > segF_dVn
const std::vector< std::vector< int > > & jacobianStamp() const
std::vector< double > segF_dVp
std::vector< double > dHQ_dH
std::vector< double > segCEquFvalue
const SolverState & solverState_
std::vector< double > segMEquQvalue
Descriptor & addPar(const char *parName, T default_value, T U::*varPtr)
Adds the parameter description to the parameter map.
Definition: N_DEV_Pars.h:1429
std::vector< double > dmF_dm
std::vector< double > daQ_da
std::vector< int > li_internalVars
std::vector< int > SegVEqnHOffset
std::vector< double > dbF_dV
std::vector< int > MEquVNodeOffset
std::vector< std::string > membraneIndpFEqus
RCP< MembraneModel > membraneModel_
static void loadModelParameters(ParametricData< Model > &model_parameters)
Linear::Vector * nextSolVectorPtr
std::vector< double > dcF_dCa
std::vector< double > segF_db
std::vector< int > CEquCNodeOffset
std::vector< int > MEquMNodeOffset
bool updateTemperature(const double &temp_tmp)
Linear::Vector * daeQVectorPtr
std::vector< double > daF_dV
Pure virtual class to augment a linear system.
void addInternalNode(Util::SymbolTable &symbol_table, int index, const InstanceName &instance_name, const std::string &lead_name)
bool processParams()
processParams
std::vector< double > segF_dm
std::vector< double > dbQ_db
std::vector< double > dcF_dV
std::vector< std::string > membraneCurrentEqus
std::vector< int > CaEquH_NodeOffset
std::vector< int > M_EquVNodeOffset
#define AssertLIDs(cmp)
std::vector< double > segHEquQvalue
std::map< int, int > segMap
virtual void registerJacLIDs(const JacobianStamp &jacLIDVec)
std::vector< double > dmQ_dm
std::vector< double > segCEquQvalue
std::vector< double > dnF_dn
std::vector< double > dCaF_dV
std::vector< double > segHEquFvalue
static Device * factory(const Configuration &configuration, const FactoryBlock &factory_block)
std::vector< std::string > membraneFunctions
std::vector< int > HEquHNodeOffset
std::vector< double > segF_dc
std::vector< int > H_EquVNodeOffset
std::vector< double > dCaQ_dCa
std::vector< double > segNEquQvalue
std::vector< int > CaEquCaNodeOffset
std::vector< double > segF_dH
std::vector< double > segCaEquQvalue
std::vector< Param > params
Parameters from the line.
void setParams(const std::vector< Param > &params)
const std::string & getName() const
static void loadInstanceParameters(ParametricData< Instance > &instance_parameters)
Definition: N_DEV_Neuron6.C:79
std::vector< int > HEquVNodeOffset
std::vector< int > CaEquM_NodeOffset
std::vector< double > segBEquQvalue
The FactoryBlock contains parameters needed by the device, instance and model creation functions...
std::vector< double > segAEquFvalue
std::vector< double > segM_EquFvalue
std::vector< int > BEquBNodeOffset
void registerJacLIDs(const std::vector< std::vector< int > > &jacLIDVec)
std::vector< double > segMEquFvalue
std::vector< std::string > membraneParameters
std::vector< double > dbF_db
std::vector< double > dnF_dV
std::vector< double > potassiumCurrent
std::vector< double > segCaEquFvalue
std::vector< int > SegVEqnVpreOffset
void registerLIDs(const std::vector< int > &intLIDVecRef, const std::vector< int > &extLIDVecRef)
std::vector< Instance * > instanceContainer
const DeviceOptions & deviceOptions_
void varTypes(std::vector< char > &varTypeVec)
std::vector< double > dhQ_dh
std::map< int, int > prevMap
std::vector< double > segQvalue
std::vector< double > dMF_dV
static Config< T > & addConfiguration()
Adds the device to the Xyce device configuration.
std::vector< double > segF_dh
std::vector< int > AEquVNodeOffset
std::vector< double > segF_dn
Linear::Matrix * dFdxMatrixPtr
std::vector< int > CEquVNodeOffset
std::vector< int > li_NaCurrentState
The Device class is an interface for device implementations.
Definition: N_DEV_Device.h:101
std::vector< std::vector< int > > jacStamp
std::vector< int > SegVEqnVsegOffset
std::vector< double > segM_EquQvalue
std::vector< double > dhF_dh
void loadNodeSymbols(Util::SymbolTable &symbol_table) const
Populates and returns the store name map.
std::vector< int > SegVEqnMOffset
std::vector< double > dhF_dV
Class Configuration contains device configuration data.
void registerStateLIDs(const std::vector< int > &staLIDVecRef)
std::vector< double > dMQ_dM
std::vector< int > H_EquH_NodeOffset
std::vector< std::vector< int > > jacobianOffsets
std::vector< double > segF_da
std::vector< double > segAEquQvalue
std::vector< double > dcQ_dc
std::vector< double > dCaF_dM
std::vector< int > CEquCaNodeOffset
const SolverState & getSolverState() const
std::vector< double > dmF_dV
#define M_PI
std::vector< int > CaEquVNodeOffset
std::vector< double > dHF_dV
std::vector< int > SegVEqnVnexOffset
std::vector< int > NEquVNodeOffset
std::vector< double > dcF_dc
virtual void forEachInstance(DeviceInstanceOp &op) const
Apply a device instance "op" to all instances associated with this model.
std::vector< double > segH_EquQvalue
std::vector< double > segH_EquFvalue
std::vector< int > li_KCurrentState
std::vector< int > BEquVNodeOffset
Linear::Vector * daeFVectorPtr
std::vector< double > dCaF_dCa
std::vector< double > segNEquFvalue
std::vector< double > dMF_dM
std::vector< double > segF_dV
std::vector< int > M_EquM_NodeOffset
std::vector< std::string > membraneIndpVars
std::vector< double > segFvalue
std::vector< double > daF_da
std::vector< std::string > membraneIndpQEqus
Instance(const Configuration &configuration, const InstanceBlock &IB, Model &Miter, const FactoryBlock &factory_block)
bool processInstanceParams()
processInstanceParams
std::vector< double > dnQ_dn
std::vector< int > AEquANodeOffset
ModelBlock represents a .MODEL line from the netlist.
Parameter may be specified as a solution dependent expression from netlist.
Definition: N_DEV_Pars.h:68
std::vector< double > sodiumCurrent
Manages parameter binding for class C.
Definition: N_DEV_Pars.h:214
InstanceBlock represent a device instance line from the netlist.
std::vector< Param > params
std::vector< double > segBEquFvalue
Linear::Matrix * dQdxMatrixPtr
void setModParams(const std::vector< Param > &params)
std::vector< double > dHF_dH
std::vector< int > NEquNNodeOffset