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