Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_DEV_Neuron3.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_Neuron3.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.33.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 #ifdef HAVE_CMATH
51 #include <cmath>
52 #else
53 #include <math.h>
54 #endif
55 
56 #include <N_UTL_Misc.h>
57 
58 // ---------- Xyce Includes ----------
59 #include <N_DEV_DeviceOptions.h>
60 #include <N_DEV_DeviceMaster.h>
61 #include <N_DEV_ExternData.h>
62 #include <N_DEV_MatrixLoadData.h>
63 #include <N_DEV_Neuron3.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 Neuron3 {
78 
79 
81 {
82 // Set up map for double precision variables:
83  p.addPar ("R", 0.0, false, ParameterType::NO_DEP,
85  &Neuron3::Instance::rIntGiven, U_OHMMM1, CAT_NONE, "Intracellular resistivity");
86 
87  p.addPar ("A", 0.0, false, ParameterType::NO_DEP,
89  &Neuron3::Instance::radiusGiven, U_METER, CAT_NONE, "Segment radius");
90 
91  p.addPar ("L", 0.0, false, ParameterType::NO_DEP,
94 
95  p.addPar ("RPS", 1.0e-6, false, ParameterType::NO_DEP,
97  &Neuron3::Instance::rIntPreviousGiven, U_OHMMM1, CAT_NONE, "Previous segment, intracellular resistivity");
98 
99  p.addPar ("APS", 0.0, false, ParameterType::NO_DEP,
101  &Neuron3::Instance::radiusPreviousGiven, U_METER, CAT_NONE, "Previous segment, segment radius");
102 
103  p.addPar ("LPS", 0.0, false, ParameterType::NO_DEP,
105  &Neuron3::Instance::lengthPreviousGiven, U_METER, CAT_NONE, "Previous segment length");
106 
107  p.addPar ("RNS", 1.0e-6, false, ParameterType::NO_DEP,
109  &Neuron3::Instance::rIntNextGiven, U_OHMMM1, CAT_NONE, "Next segment, intracellular resistivity");
110 
111  p.addPar ("ANS", 0.0, false, ParameterType::NO_DEP,
113  &Neuron3::Instance::radiusNextGiven, U_METER, CAT_NONE, "Next segment, segment radius");
114 
115  p.addPar ("LNS", 0.0, false, ParameterType::NO_DEP,
117  &Neuron3::Instance::lengthNextGiven, U_METER, CAT_NONE, "Next segment length");
118 
119  // add non-double types
121  &Neuron3::Instance::nSegGiven , U_NONE, CAT_NONE, "Number of segments" );
122 }
123 
125 {
126  // Set up map for normal (double) param variables:
127  p.addPar ("CMEM", 0.0, false, ParameterType::NO_DEP,
130  U_FARAD, CAT_NONE, "Membrane capacitance");
131 
132  p.addPar ("GMEM", 0.0, false, ParameterType::NO_DEP,
135  U_OHMM1, CAT_NONE, "Membrane conductance");
136 
137  p.addPar ("VREST", 0.0, false, ParameterType::NO_DEP,
140  U_VOLT, CAT_NONE, "Resting potential");
141 
142  p.addPar ("EK", 0.0, false, ParameterType::NO_DEP,
145  U_VOLT, CAT_NONE, "Potassium resting potential");
146 
147  p.addPar ("GK", 0.0, false, ParameterType::NO_DEP,
150  U_OHMM1, CAT_NONE, "Potassium base conductance");
151 
152  p.addPar ("ENA", 0.0, false, ParameterType::NO_DEP,
155  U_VOLT, CAT_NONE, "Sodium resting potential");
156 
157  p.addPar ("GNA", 0.0, false, ParameterType::NO_DEP,
160  U_OHMM1, CAT_NONE, "Sodium base conductance");
161 
162  p.addPar ("R", 0.0, false, ParameterType::NO_DEP,
165  U_OHMMM1, CAT_NONE, "Intracellular resistivity");
166 
167  p.addPar ("A", 0.0, false, ParameterType::NO_DEP,
170  U_METER, CAT_NONE, "Segment radius");
171 
172  p.addPar ("L", 0.0, false, ParameterType::NO_DEP,
175  U_METER, CAT_NONE, "Cable length");
176 
177  p.addPar ("RPS", 1.0e-6, false, ParameterType::NO_DEP,
180  U_OHMMM1, CAT_NONE, "Previous segment, intracellular resistivity");
181 
182  p.addPar ("APS", 0.0, false, ParameterType::NO_DEP,
185  U_METER, CAT_NONE, "Previous segment, segment radius");
186 
187  p.addPar ("LPS", 0.0, false, ParameterType::NO_DEP,
190  U_METER, CAT_NONE, "Previous segment length");
191 
192  p.addPar ("RNS", 1.0e-6, false, ParameterType::NO_DEP,
195  U_OHMMM1, CAT_NONE, "Next segment, intracellular resistivity");
196 
197  p.addPar ("ANS", 0.0, false, ParameterType::NO_DEP,
200  U_METER, CAT_NONE, "Next segment, segment radius");
201 
202  p.addPar ("LNS", 0.0, false, ParameterType::NO_DEP,
205  U_METER, CAT_NONE, "Next segment length");
206 
207  // add non-double types
209  &Neuron3::Model::nSegGiven , U_NONE, CAT_NONE, "Number of segments" );
210 }
211 
212 
213 
214 //
215 // static class member inits
216 //
217 
218 
219 // Class Instance
220 //-----------------------------------------------------------------------------
221 // Function : Instance::processParams
222 // Purpose :
223 // Special Notes :
224 // Scope : public
225 // Creator : Richard Schiek, Electrical and Microsytem Modeling
226 // Creation Date : 06/10/09
227 //-----------------------------------------------------------------------------
229 {
230  // If there are any time dependent parameters, set their values at for
231  // the current time.
232 
233  // now set the temperature related stuff.
234  //updateTemperature(temp);
235 
236  return true;
237 }
238 
239 //-----------------------------------------------------------------------------
240 // Function : Instance::updateTemperature
241 // Purpose :
242 // Special Notes :
243 // Scope : public
244 // Creator : Richard Schiek, Electrical and Microsytem Modeling
245 // Creation Date : 06/10/09
246 //-----------------------------------------------------------------------------
247 bool Instance::updateTemperature ( const double & temp)
248 {
249  bool bsuccess = true;
250  return bsuccess;
251 }
252 
253 //-----------------------------------------------------------------------------
254 // Function : Instance::Instance
255 // Purpose : constructor
256 // Special Notes :
257 // Scope : public
258 // Creator : Richard Schiek, Electrical and Microsytem Modeling
259 // Creation Date : 06/10/09
260 //-----------------------------------------------------------------------------
262  const Configuration & configuration,
263  const InstanceBlock & IB,
264  Model & Miter,
265  const FactoryBlock & factory_block)
266  : DeviceInstance(IB, configuration.getInstanceParameters(), factory_block),
267  model_(Miter),
268  rInt(0.0),
269  radius(0.0),
270  length(0.0),
271  segArea(0.0),
272  nSeg(0),
273  rIntGiven(false),
274  radiusGiven(false),
275  lengthGiven(false),
276  nSegGiven(false),
277  rIntPrevious(0.0),
278  radiusPrevious(0.0),
279  lengthPrevious(0.0),
280  rIntNext(0.0),
281  radiusNext(0.0),
282  lengthNext(0.0),
283  rIntPreviousGiven(false),
284  radiusPreviousGiven(false),
285  lengthPreviousGiven(false),
286  rIntNextGiven(false),
287  radiusNextGiven(false),
288  lengthNextGiven(false),
289  kcl1Fvalue(0.0),
290  kcl2Fvalue(0.0),
291  dkcl1F_dVin(0.0),
292  dkcl1F_dVs0(0.0),
293  dkcl2F_dVout(0.0),
294  dkcl2F_dVsn(0.0),
295  li_Pos(0),
296  li_Neg(0),
297  APosEquPosNodeOffset(0),
298  APosEquNextNodeOffset(0),
299  ANegEquNegNodeOffset(0),
300  ANegEquLastNodeOffset(0)
301 {
302  // we have pased the model and instance parameters so now we can calculate
303  // the number of internal vars
304  numExtVars = 2; // input and output voltage
305 
306  // Set params to constant default values:
307  setDefaultParams ();
308 
309  // Set params according to instance line
310  setParams (IB.params);
311 
312  // Set any non-constant parameter defaults:
313 
314  //if (!given("TEMP"))
315  // temp = getDeviceOptions().temp.dVal();
316 
317  // Calculate any parameters specified as expressions:
319 
320  // calculate dependent (ie computed) params and check for errors:
321  processParams ();
322 
323  // pull unspecified params out of the model if they weren't specified here
324  if( !nSegGiven && model_.nSegGiven )
325  {
326  nSeg = model_.nSeg;
327  nSegGiven = true;
328  }
329  if( !rIntGiven && model_.rIntGiven )
330  {
331  rInt = model_.rInt;
332  rIntGiven = true;
333  }
334  if( !radiusGiven && model_.radiusGiven )
335  {
336  radius = model_.radius;
337  radiusGiven = true;
338  }
339  if( !lengthGiven && model_.lengthGiven )
340  {
341  length = model_.length;
342  lengthGiven = true;
343  }
344 
345  // if nSeg is still unknown then estimate it via lambda-d rule (TO DO)
346  if( !nSegGiven )
347  {
348  nSeg = 10;
349  }
350 
351  // now that nSeg, length and radius are set calculate segment area
352  segArea = 2.0 * M_PI * radius * length / nSeg;
353 
354  // now we can calculate the number of internal vars
355  numIntVars = nSeg * 4; // ion channel vars + one internal voltage node var for each segment
356  numStateVars = nSeg * 2;
357  int numVars = numExtVars + numIntVars;
358 
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
362  //
363  // Vin : g(0,1) * (V(1) - Vin ) = 0
364  // Vout : g(n,n-1) * (V(n-1) - Vout) = 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
369  // Vin Vout V1 n m h V2 n m h V(nSeg) n m h
370  // kcl Vin -g(0,1) g(0,1)
371  // kcl Vout -g(n,n-1) g(n,n-1)
372  // kcl V1 yes yes yes yes yes yes
373  // a yes yes
374  // b yes yes
375  // m yes yes
376  // kcl V2 yes yes yes yes yes yes
377  // a 2 yes yes
378  // b 2 yes yes
379  // m 2 yes yes
380  // kcl VnSeg yes yes yes yes yes yes
381  // a nSeg yes yes
382  // b nSeg yes yes
383  // m nSeg yes yes
384  //
385  // jacobian element count by row:
386  // vin: 2
387  // vout: 2
388  // v1: 6
389  // a1: 2
390  // b1: 2
391  // m1: 2
392  // v2: 6
393  // a2: 2
394  // b2: 2
395  // m2: 2
396  // vnSeg:6
397  // a: 2
398  // b: 2
399  // m: 2
400  //
401  // set up jacStamp
402  if( jacStamp.empty() ) // redundant as jacStamp is not static for this device
403  { // it can't be as each cable may have a different number of nodes
404  jacStamp.resize(numVars);
405  jacStamp[0].resize(2);
406  jacStamp[0][0] = 0;
407  jacStamp[0][1] = 1;
408  jacStamp[1].resize(2);
409  jacStamp[1][0] = 1;
410  jacStamp[1][1] = numVars - 4;
411  for( int i=2; i<numVars; i+=4)
412  {
413  jacStamp[i].resize(6);
414  if( i == 2 )
415  {
416  jacStamp[i][0] = 0;
417  }
418  else
419  {
420  jacStamp[i][0] = i-4;
421  }
422  jacStamp[i][1] = i;
423  jacStamp[i][2] = i+1;
424  jacStamp[i][3] = i+2;
425  jacStamp[i][4] = i+3;
426  if( i==(numVars-4) )
427  {
428  jacStamp[i][5] = 1;
429  }
430  else
431  {
432  jacStamp[i][5] = i+4;
433  }
434  jacStamp[i+1].resize(2);
435  jacStamp[i+1][0] = i;
436  jacStamp[i+1][1] = i+1;
437  jacStamp[i+2].resize(2);
438  jacStamp[i+2][0] = i;
439  jacStamp[i+2][1] = i+2;
440  jacStamp[i+3].resize(2);
441  jacStamp[i+3][0] = i;
442  jacStamp[i+3][1] = i+3;
443  }
444 
445  }
446 
447  /*
448  // print out jacStamp
449  int numRows = jacStamp.size();
450  for( int i=0; i< numRows; i++ )
451  {
452  int numCol = jacStamp[i].size();
453  Xyce::dout() << "jacStamp[ " << i << " ] = { ";
454  for(int j=0; j<numCol; j++)
455  {
456  Xyce::dout() << jacStamp[i][j] << " ";
457  }
458  Xyce::dout() << " } " << std::endl;
459  }
460  Xyce::dout() << std::endl;
461  */
462 
463  // calculate segment conductivities used in load calls:
464  // Note: conductivity to the previous and next segment is not symmetric if the segment length and/or radius is not are not equal
465  // the full formula is:
466  // g(n,n') = radius * (radius')^2 / ( rInt segLength * ( segLength * (radius')^2 + segLength' * (radius)^2 ) )
467  //
468  // equation 6.30 Theoretical neuroscience: computational and mathematical modeling of neural systems, Peter Dayan and L. F. Abbot 2001
469  // If we allow variable segment lengths and radii then we'll need to update this calculation
470  //
471  // Note: the above equation has the wrong units unless rInt (which is Ohm Length) is rLong (which is Ohm/Length).
472  //
473  gForward.resize(nSeg);
474  gBackward.resize(nSeg);
475  double segLength = length / nSeg;
476  double rLong = rInt / (M_PI * radius * radius); // longitudinal resistivity (ohm/length)
477  gBackward[0] = radius * (radiusPrevious * radiusPrevious) / (rLong * segLength * ( segLength * radiusPrevious * radiusPrevious + lengthPrevious * radius * radius ));
478  gForward[0] = radius * (radius * radius) / (rLong * segLength * ( segLength * radius * radius + segLength * radius * radius ));
479  gBackward[nSeg-1] = radius * (radius * radius) / (rLong * segLength * ( segLength * radius * radius + segLength * radius * radius ));
480  gForward[nSeg-1] = radius * (radiusNext * radiusNext) / (rLong * segLength * ( segLength * radiusNext * radiusNext + lengthNext * radius * radius ));
481  for(int i=1; i<(nSeg-1); i++)
482  {
483  gBackward[i] = radius * (radius * radius) / (rLong * segLength * ( segLength * radius * radius + segLength * radius * radius ));
484  gForward[i] = gBackward[i];
485  }
486  /*
487  gBackward[0] = radius * (radiusPrevious * radiusPrevious) / (rInt * segLength * ( segLength * radiusPrevious * radiusPrevious + lengthPrevious * radius * radius ));
488  gForward[0] = radius * (radius * radius) / (rInt * segLength * ( segLength * radius * radius + segLength * radius * radius ));
489  gBackward[nSeg-1] = radius * (radius * radius) / (rInt * segLength * ( segLength * radius * radius + segLength * radius * radius ));
490  gForward[nSeg-1] = radius * (radiusNext * radiusNext) / (rInt * segLength * ( segLength * radiusNext * radiusNext + lengthNext * radius * radius ));
491  for(int i=1; i<(nSeg-1); i++)
492  {
493  gBackward[i] = radius * (radius * radius) / (rInt * segLength * ( segLength * radius * radius + segLength * radius * radius ));
494  gForward[i] = gBackward[i];
495  }
496  */
497  // allocate space for load and jacobian terms per segment
498  // variable indecies loads
499  li_Vol.resize(nSeg);
500  li_nPro.resize(nSeg);
501  li_mPro.resize(nSeg);
502  li_hPro.resize(nSeg);
503  li_KCurrentState.resize(nSeg);
504  li_NaCurrentState.resize(nSeg);
505  segFvalue.resize(nSeg);
506  segQvalue.resize(nSeg);
507  segNEquFvalue.resize(nSeg);
508  segNEquQvalue.resize(nSeg);
509  segMEquFvalue.resize(nSeg);
510  segMEquQvalue.resize(nSeg);
511  segHEquFvalue.resize(nSeg);
512  segHEquQvalue.resize(nSeg);
513 
514  // jacobian elements
515  segF_dVp.resize(nSeg);
516  segF_dV.resize(nSeg);
517  segF_dVn.resize(nSeg);
518  segF_dn.resize(nSeg);
519  segF_dm.resize(nSeg);
520  segF_dh.resize(nSeg);
521  segQ_dV.resize(nSeg);
522  dnF_dV.resize(nSeg);
523  dnF_dn.resize(nSeg);
524  dnQ_dn.resize(nSeg);
525  dmF_dV.resize(nSeg);
526  dmF_dm.resize(nSeg);
527  dmQ_dm.resize(nSeg);
528  dhF_dV.resize(nSeg);
529  dhF_dh.resize(nSeg);
530  dhQ_dh.resize(nSeg);
531 
532  // state vars
533  potassiumCurrent.resize(nSeg);
534  sodiumCurrent.resize(nSeg);
535 
536 }
537 
538 //-----------------------------------------------------------------------------
539 // Function : Instance::~Instance
540 // Purpose : destructor
541 // Special Notes :
542 // Scope : public
543 // Creator : Richard Schiek, Electrical and Microsytem Modeling
544 // Creation Date : 06/10/09
545 //-----------------------------------------------------------------------------
547 {
548 }
549 
550 //-----------------------------------------------------------------------------
551 // Function : Instance::registerLIDs
552 // Purpose :
553 // Special Notes :
554 // Scope : public
555 // Creator : Richard Schiek, Electrical and Microsytem Modeling
556 // Creation Date : 06/10/09
557 //-----------------------------------------------------------------------------
558 void Instance::registerLIDs(const std::vector<int> & intLIDVecRef,
559  const std::vector<int> & extLIDVecRef)
560 {
561  AssertLIDs(intLIDVecRef.size() == numIntVars);
562  AssertLIDs(extLIDVecRef.size() == numExtVars);
563 
564 #ifdef Xyce_DEBUG_DEVICE
565  if (getDeviceOptions().debugLevel > 0)
566  {
567  Xyce::dout() << std::endl << section_divider << std::endl;
568  Xyce::dout() << " Instance::registerLIDs" << std::endl;
569  Xyce::dout() << " name = " << getName() << std::endl;
570  }
571 #endif
572 
573  // copy over the global ID lists.
574  intLIDVec = intLIDVecRef;
575  extLIDVec = extLIDVecRef;
576 
577  li_Pos = extLIDVec[0];
578  li_Neg = extLIDVec[1];
579  for( int i=0, j=0; i<nSeg; i++, j+=4)
580  {
581  li_Vol[i] = intLIDVec[j];
582  li_nPro[i] = intLIDVec[j+1];
583  li_mPro[i] = intLIDVec[j+2];
584  li_hPro[i] = intLIDVec[j+3];
585  }
586 
587 #ifdef Xyce_DEBUG_DEVICE
588  if (getDeviceOptions().debugLevel > 0 )
589  {
590  Xyce::dout() << " li_Pos = " << li_Pos << std::endl
591  << " li_Neg = " << li_Neg << std::endl;
592  for( int i=0; i<nSeg; i++ )
593  {
594  Xyce::dout() << " li_Vol[ " << i << " ] = " << li_Vol[i] << std::endl
595  << " li_nPro[ " << i << " ] = " << li_nPro[i] << std::endl
596  << " li_mPro[ " << i << " ] = " << li_mPro[i] << std::endl
597  << " li_hPro[ " << i << " ] = " << li_hPro[i] << std::endl;
598  }
599  }
600 #endif
601 
602 #ifdef Xyce_DEBUG_DEVICE
603  if (getDeviceOptions().debugLevel > 0 )
604  {
605  Xyce::dout() << section_divider << std::endl;
606  }
607 #endif
608 }
609 
610 //-----------------------------------------------------------------------------
611 // Function : Instance::getIntNameMap
612 // Purpose :
613 // Special Notes :
614 // Scope : public
615 // Creator : Richard Schiek, Electrical and Microsytem Modeling
616 // Creation Date : 06/10/09
617 //-----------------------------------------------------------------------------
618 std::map<int,std::string> & Instance::getIntNameMap ()
619 {
620  // set up the internal name map, if it hasn't been already.
621  if (intNameMap.empty ())
622  {
623  std::string tmpstr;
624  for( int i=0; i<nSeg; i++)
625  {
626  std::ostringstream segNumber;
627  segNumber << i;
628 
629  tmpstr = getName() + "_" + "V" + segNumber.str();
630  spiceInternalName (tmpstr);
631  intNameMap[ li_Vol[i] ] = tmpstr;
632 
633  tmpstr = getName() + "_" + "N" + segNumber.str();
634  spiceInternalName (tmpstr);
635  intNameMap[ li_nPro[i] ] = tmpstr;
636 
637  tmpstr = getName() + "_" + "M" + segNumber.str();
638  spiceInternalName (tmpstr);
639  intNameMap[ li_mPro[i] ] = tmpstr;
640 
641  tmpstr = getName() + "_" + "H" + segNumber.str();
642  spiceInternalName (tmpstr);
643  intNameMap[ li_hPro[i] ] = tmpstr;
644  }
645  }
646  return intNameMap;
647 }
648 
649 //-----------------------------------------------------------------------------
650 // Function : Instance::registerStateLIDs
651 // Purpose :
652 // Special Notes :
653 // Scope : public
654 // Creator : Richard Schiek, Electrical and Microsytem Modeling
655 // Creation Date : 06/10/09
656 //-----------------------------------------------------------------------------
657 void Instance::registerStateLIDs( const std::vector<int> & staLIDVecRef )
658 {
659  AssertLIDs(staLIDVecRef.size() == numStateVars);
660 
661  // copy over the global ID lists.
662  staLIDVec = staLIDVecRef;
663 
664  for( int i=0, j=0; i<nSeg; i++, j+=2)
665  {
666  li_KCurrentState[i] = staLIDVec[j];
667  li_NaCurrentState[i] = staLIDVec[j+1];
668  }
669 
670 }
671 
672 //-----------------------------------------------------------------------------
673 // Function : Instance::loadDeviceMask
674 //
675 // Purpose : Loads the zero elements of the device mask
676 //
677 // Special Notes : elements of the error vector associated with zero
678 // elements of the mask will not be included in weighted
679 // norms by the time integrator.
680 //
681 // Scope : public
682 // Creator : Tom Russo, SNL, Electrical and Microsystems Modeling
683 // Creation Date : 06/10/09
684 //-----------------------------------------------------------------------------
686 {
687  bool returnVal=false;
688  N_LAS_Vector * maskVectorPtr = extData.deviceMaskVectorPtr;
689 
690 // Xyce::dout() << "Masking n, m and h" << std::endl;
691 // (*maskVectorPtr)[li_nPro] = 0.0;
692 // (*maskVectorPtr)[li_mPro] = 0.0;
693 // (*maskVectorPtr)[li_hPro] = 0.0;
694 // returnVal = true;
695 
696  return (returnVal);
697 }
698 
699 //-----------------------------------------------------------------------------
700 // Function : Instance::jacobianStamp
701 // Purpose :
702 // Special Notes :
703 // Scope : public
704 // Creator : Richard Schiek, Electrical and Microsytem Modeling
705 // Creation Date : 06/10/09
706 //-----------------------------------------------------------------------------
707 const std::vector< std::vector<int> > & Instance::jacobianStamp() const
708 {
709  return jacStamp;
710 }
711 
712 //-----------------------------------------------------------------------------
713 // Function : Instance::registerJacLIDs
714 // Purpose :
715 // Special Notes :
716 // Scope : public
717 // Creator : Richard Schiek, Electrical and Microsytem Modeling
718 // Creation Date : 06/10/09
719 //-----------------------------------------------------------------------------
720 void Instance::registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec )
721 {
722  DeviceInstance::registerJacLIDs( jacLIDVec );
723 
724  // external terminals
725  APosEquPosNodeOffset = jacLIDVec[0][0];
726  APosEquNextNodeOffset = jacLIDVec[0][1];
727  ANegEquNegNodeOffset = jacLIDVec[1][0];
728  ANegEquLastNodeOffset = jacLIDVec[1][1];
729  /*
730  Xyce::dout() << "APosEquPosNodeOffset = " << APosEquPosNodeOffset << std::endl;
731  Xyce::dout() << "APosEquNextNodeOffset = " << APosEquNextNodeOffset << std::endl;
732  Xyce::dout() << "ANegEquNegNodeOffset = " << ANegEquNegNodeOffset << std::endl;
733  Xyce::dout() << "ANegEquLastNodeOffset = " << ANegEquLastNodeOffset << std::endl;
734  */
735 
736  // internal variables
737  SegVEqnVpreOffset.resize(nSeg);
738  SegVEqnVsegOffset.resize(nSeg);
739  SegVEqnVnexOffset.resize(nSeg);
740  SegVEqnNOffset.resize(nSeg);
741  SegVEqnMOffset.resize(nSeg);
742  SegVEqnHOffset.resize(nSeg);
743  NEquVNodeOffset.resize(nSeg);
744  NEquNNodeOffset.resize(nSeg);
745  MEquVNodeOffset.resize(nSeg);
746  MEquMNodeOffset.resize(nSeg);
747  HEquVNodeOffset.resize(nSeg);
748  HEquHNodeOffset.resize(nSeg);
749  for(int i=0, j=2; i<nSeg; i++, j+=4 )
750  {
751  SegVEqnVpreOffset[i] = jacLIDVec[j][0];
752  SegVEqnVsegOffset[i] = jacLIDVec[j][1];
753  SegVEqnNOffset[i] = jacLIDVec[j][2];
754  SegVEqnMOffset[i] = jacLIDVec[j][3];
755  SegVEqnHOffset[i] = jacLIDVec[j][4];
756  SegVEqnVnexOffset[i] = jacLIDVec[j][5];
757 
758  NEquVNodeOffset[i] = jacLIDVec[j+1][0];
759  NEquNNodeOffset[i] = jacLIDVec[j+1][1];
760  MEquVNodeOffset[i] = jacLIDVec[j+2][0];
761  MEquMNodeOffset[i] = jacLIDVec[j+2][1];
762  HEquVNodeOffset[i] = jacLIDVec[j+3][0];
763  HEquHNodeOffset[i] = jacLIDVec[j+3][1];
764  /*
765  Xyce::dout() << SegVEqnVpreOffset[i] << ", "
766  << SegVEqnVsegOffset[i] << ", "
767  << SegVEqnNOffset[i] << ", "
768  << SegVEqnMOffset[i] << ", "
769  << SegVEqnHOffset[i] << ", "
770  << SegVEqnVnexOffset[i] << ", "
771  << NEquVNodeOffset[i] << ", "
772  << NEquNNodeOffset[i] << ", "
773  << MEquVNodeOffset[i] << ", "
774  << MEquMNodeOffset[i] << ", "
775  << HEquVNodeOffset[i] << ", "
776  << HEquHNodeOffset[i] << std::endl;
777  */
778  }
779 }
780 
781 
782 //-----------------------------------------------------------------------------
783 // Function : Instance::updateIntermediateVars
784 // Purpose :
785 // Special Notes :
786 // Scope : public
787 // Creator : Richard Schiek, Electrical and Microsytem Modeling
788 // Creation Date : 06/10/09
789 //-----------------------------------------------------------------------------
791 {
792  //Xyce::dout() << "Instance::updateIntermediateVars()" << std::endl;
793 
794  bool bsuccess = true;
795 
796  // here we take the current solutions for V1, V2, n, m and h
797  // and use those to calculate all the terms needed for the next
798  // load cycle (F, Q, dFdX, dQdX)
799  N_LAS_Vector * solVectorPtr = extData.nextSolVectorPtr;
800 
801  double vIn = (*solVectorPtr)[li_Pos];
802  double vOut = (*solVectorPtr)[li_Neg];
803 
804  // take care of the input and output nodes as they are different
805  kcl1Fvalue = gForward[0] * ((*solVectorPtr)[li_Vol[0]] - vIn );
806  dkcl1F_dVin = -gForward[0];
807  dkcl1F_dVs0 = gForward[0];
808  kcl2Fvalue = gBackward[nSeg-1] * ((*solVectorPtr)[li_Vol[nSeg-1]] - vOut );
811  //Xyce::dout() << "end update: " << vIn << ", " << vOut << ", " << gForward[0] << ", " << gBackward[nSeg-1] << std::endl;
812 
813  // loop over segments getting all the load and jacobian terms for each segment
814  for( int i=0; i<nSeg; i++ )
815  {
816  // for this segment get the values of the local vars
817  double vSeg = (*solVectorPtr)[li_Vol[i]];
818  double vNext = 0.0;
819  if (i == (nSeg - 1))
820  {
821  vNext = vOut;
822  }
823  else
824  {
825  vNext = (*solVectorPtr)[li_Vol[i+1]];
826  }
827  double vPrev = 0.0;
828  if (i == 0 )
829  {
830  vPrev = vIn;
831  }
832  else
833  {
834  vPrev = (*solVectorPtr)[li_Vol[i-1]];
835  }
836  double nVarSeg = (*solVectorPtr)[li_nPro[i]];
837  double mVarSeg = (*solVectorPtr)[li_mPro[i]];
838  double hVarSeg = (*solVectorPtr)[li_hPro[i]];
839 
840  // do the voltage equation for this node
841  // get function and derivative values as we go.
842  {
843  // F part
844  // use scoping to avoid lots of similar variable names
845  const int numDeriv = 6;
846  Sacado::Fad::SFad<double,6> vVar( numDeriv, 0, vSeg );
847  Sacado::Fad::SFad<double,6> vVpr( numDeriv, 1, vPrev );
848  Sacado::Fad::SFad<double,6> vVne( numDeriv, 2, vNext );
849  Sacado::Fad::SFad<double,6> nVar( numDeriv, 3, nVarSeg );
850  Sacado::Fad::SFad<double,6> mVar( numDeriv, 4, mVarSeg );
851  Sacado::Fad::SFad<double,6> hVar( numDeriv, 5, hVarSeg );
852  // parameters
853  Sacado::Fad::SFad<double,6> gPrev( gBackward[i] );
854  Sacado::Fad::SFad<double,6> gNext( gForward[i] );
855  Sacado::Fad::SFad<double,6> gMemVar( model_.gMem * segArea );
856  Sacado::Fad::SFad<double,6> vRestVar( model_.vRest );
857  Sacado::Fad::SFad<double,6> gKVar( model_.gK * segArea );
858  Sacado::Fad::SFad<double,6> eKVar( model_.eK );
859  Sacado::Fad::SFad<double,6> gNaVar( model_.gNa * segArea );
860  Sacado::Fad::SFad<double,6> eNaVar( model_.eNa );
861  //Xyce::dout() << "segment update: " << i << ", " << vSeg << ", " << vPrev << ", " << vNext << ", " << gBackward[i] << ", " << gForward[i] << std::endl;
862  // compute the value and derivative terms for KCL 1 F
863  Sacado::Fad::SFad<double,6> resultFad;
864  resultFad = kcl1EquF( vVar, vVpr, vVne, nVar, mVar, hVar, gPrev, gNext, gMemVar, vRestVar, gKVar, eKVar, gNaVar, eNaVar );
865  segFvalue[i] = resultFad.val();
866  segF_dV[i] = resultFad.dx(0);
867  segF_dVp[i] = resultFad.dx(1);
868  segF_dVn[i] = resultFad.dx(2);
869  segF_dn[i] = resultFad.dx(3);
870  segF_dm[i] = resultFad.dx(4);
871  segF_dh[i] = resultFad.dx(5);
872  }
873  {
874  // Q part
875  const int numDeriv = 2;
876  Sacado::Fad::SFad<double,2> vVar( numDeriv, 0, vSeg );
877 
878  // parameters
879  Sacado::Fad::SFad<double,2> cMemVar( model_.cMem * segArea );
880 
881  Sacado::Fad::SFad<double,2> resultFad;
882  resultFad = kcl1EquQ( vVar, cMemVar );
883  segQvalue[i] = resultFad.val();
884  segQ_dV[i] = resultFad.dx(0);
885 
886  }
887 
888  // n - equation
889  {
890  const int numDeriv = 2;
891  Sacado::Fad::SFad<double,2> vVar( numDeriv, 0, vSeg );
892  Sacado::Fad::SFad<double,2> nVar( numDeriv, 1, nVarSeg );
893  // parameter
894  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
895 
896  Sacado::Fad::SFad<double,2> resultFad = nEquF( vVar, nVar, vRestVar);
897  segNEquFvalue[i] = resultFad.val();
898  dnF_dV[i] = resultFad.dx(0);
899  dnF_dn[i] = resultFad.dx(1);
900  }
901  {
902  const int numDeriv = 1;
903  Sacado::Fad::SFad<double,1> nVar( numDeriv, 0, nVarSeg );
904 
905  Sacado::Fad::SFad<double,1> resultFad = nEquQ( nVar );
906  segNEquQvalue[i] = resultFad.val();
907  dnQ_dn[i] = resultFad.dx(0);
908  }
909 
910  // m - equation
911  {
912  const int numDeriv = 2;
913  Sacado::Fad::SFad<double,2> vVar( numDeriv, 0, vSeg );
914  Sacado::Fad::SFad<double,2> mVar( numDeriv, 1, mVarSeg );
915  // parameter
916  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
917 
918  Sacado::Fad::SFad<double,2> resultFad = mEquF( vVar, mVar, vRestVar );
919  segMEquFvalue[i] = resultFad.val();
920  dmF_dV[i] = resultFad.dx(0);
921  dmF_dm[i] = resultFad.dx(1);
922  }
923  {
924  const int numDeriv = 1;
925  Sacado::Fad::SFad<double,1> mVar( numDeriv, 0, mVarSeg );
926 
927  Sacado::Fad::SFad<double,1> resultFad = mEquQ( mVar );
928  segMEquQvalue[i] = resultFad.val();
929  dmQ_dm[i] = resultFad.dx(0);
930  }
931 
932  // h - equation
933  {
934  const int numDeriv = 2;
935  Sacado::Fad::SFad<double,2> vVar( numDeriv, 0, vSeg );
936  Sacado::Fad::SFad<double,2> hVar( numDeriv, 1, hVarSeg );
937  // parameter
938  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
939 
940  Sacado::Fad::SFad<double,2> resultFad = hEquF( vVar, hVar, vRestVar );
941  segHEquFvalue[i] = resultFad.val();
942  dhF_dV[i] = resultFad.dx(0);
943  dhF_dh[i] = resultFad.dx(1);
944  }
945  {
946  const int numDeriv = 1;
947  Sacado::Fad::SFad<double,1> hVar( numDeriv, 0, hVarSeg );
948 
949  Sacado::Fad::SFad<double,1> resultFad = hEquQ( hVar );
950  segHEquQvalue[i] = resultFad.val();
951  dhQ_dh[i] = resultFad.dx(0);
952  }
953 
954  // calculate potassium current
955  potassiumCurrent[i] = model_.gK * pow(nVarSeg, 4.0) * (vSeg - model_.eK);
956 
957  // calculate sodium current
958  sodiumCurrent[i] = model_.gNa * pow(mVarSeg, 3.0) * hVarSeg * (vSeg - model_.eNa);
959 
960 #ifdef Xyce_DEBUG_DEVICE
962  {
963  Xyce::dout() << "Segment " << i << std::endl
964  << "vPrev = " << vPrev << std::endl
965  << "vSeg = " << vSeg << std::endl
966  << "vNext = " << vNext << std::endl
967  << "n, m, h = " << nVarSeg << ", " << mVarSeg << ", " << hVarSeg << std::endl
968  << "segFvalue = " << segFvalue[i] << std::endl
969  << "segQvalue = " << segQvalue[i] << std::endl
970  << "segNEquFvalue = " << segNEquFvalue[i] << std::endl
971  << "segNEquQvalue = " << segNEquQvalue[i] << std::endl
972  << "segMEquFvalue = " << segMEquFvalue[i] << std::endl
973  << "segMEquQvalue = " << segMEquQvalue[i] << std::endl
974  << "segHEquFvalue = " << segHEquFvalue[i] << std::endl
975  << "segHEquQvalue = " << segHEquQvalue[i] << std::endl
976  << std::endl;
977 
978  }
979 #endif
980 
981 
982  }
983 
984 
985 
986  return bsuccess;
987 }
988 //-----------------------------------------------------------------------------
989 // Function : Instance::updatePrimaryState
990 // Purpose :
991 // Special Notes :
992 // Scope : public
993 // Creator : Richard Schiek, Electrical and Microsytem Modeling
994 // Creation Date : 06/10/09
995 //-----------------------------------------------------------------------------
997 {
998  bool bsuccess = true;
999 
1001 
1002  N_LAS_Vector * solVectorPtr = extData.nextSolVectorPtr;
1003  N_LAS_Vector * staVectorPtr = extData.nextStaVectorPtr;
1004 
1005  for( int i=0; i<nSeg; i++)
1006  {
1007  (*staVectorPtr)[li_KCurrentState[i]] = potassiumCurrent[i];
1008  (*staVectorPtr)[li_NaCurrentState[i]] = sodiumCurrent[i];
1009  }
1010 
1011  return bsuccess;
1012 }
1013 
1014 //-----------------------------------------------------------------------------
1015 // Function : Instance::updateSecondaryState
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  N_LAS_Vector * staVectorPtr = extData.nextStaVectorPtr;
1026 
1027  return bsuccess;
1028 }
1029 
1030 //-----------------------------------------------------------------------------
1031 // Function : Instance::loadDAEQVector
1032 //
1033 // Purpose : Loads the Q-vector contributions
1034 //
1035 // Special Notes : The "Q" vector is part of a standard DAE formalism in
1036 // which the system of equations is represented as:
1037 //
1038 // f(x) = dQ(x)/dt + F(x) - B(t) = 0
1039 //
1040 // Scope : public
1041 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1042 // Creation Date : 06/10/09
1043 //-----------------------------------------------------------------------------
1045 {
1046  bool bsuccess = true;
1047 
1048  N_LAS_Vector * daeQVecPtr = extData.daeQVectorPtr;
1049 
1050  for( int i=0; i<nSeg ; i++)
1051  {
1052  (*daeQVecPtr)[li_Vol[i]] += segQvalue[i];
1053  (*daeQVecPtr)[li_nPro[i]] += segNEquQvalue[i];
1054  (*daeQVecPtr)[li_mPro[i]] += segMEquQvalue[i];
1055  (*daeQVecPtr)[li_hPro[i]] += segHEquQvalue[i];
1056  }
1057 
1058  return bsuccess;
1059 }
1060 
1061 
1062 
1063 //-----------------------------------------------------------------------------
1064 // Function : Instance::loadDAEFVector
1065 //
1066 // Purpose : Loads the F-vector contributions for a single
1067 // Neuron 3 instance.
1068 //
1069 // Special Notes :
1070 //
1071 // Scope : public
1072 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1073 // Creation Date : 06/10/09
1074 //-----------------------------------------------------------------------------
1076 {
1077  bool bsuccess=true;
1078 
1079  N_LAS_Vector * daeFVecPtr = extData.daeFVectorPtr;
1080 
1081  (*daeFVecPtr)[li_Pos] += kcl1Fvalue;
1082  (*daeFVecPtr)[li_Neg] += kcl2Fvalue;
1083 
1084  for( int i=0; i<nSeg ; i++)
1085  {
1086  (*daeFVecPtr)[li_Vol[i]] += segFvalue[i];
1087  (*daeFVecPtr)[li_nPro[i]] += segNEquFvalue[i];
1088  (*daeFVecPtr)[li_mPro[i]] += segMEquFvalue[i];
1089  (*daeFVecPtr)[li_hPro[i]] += segHEquFvalue[i];
1090  }
1091 
1092  return bsuccess;
1093 }
1094 
1095 //-----------------------------------------------------------------------------
1096 // Function : Instance::loadDAEdQdx
1097 //
1098 // Purpose : Loads the Q-vector contributions for a single
1099 // Neuron 3 instance.
1100 // Scope : public
1101 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1102 // Creation Date : 06/10/09
1103 //-----------------------------------------------------------------------------
1105 {
1106  bool bsuccess = true;
1107 
1108  N_LAS_Matrix * dQdxMatPtr = extData.dQdxMatrixPtr;
1109 
1110  for( int i=0; i<nSeg ; i++)
1111  {
1112  (*dQdxMatPtr)[li_Vol[i]][SegVEqnVsegOffset[i]] += segQ_dV[i];
1113  (*dQdxMatPtr)[li_nPro[i]][NEquNNodeOffset[i]] += dnQ_dn[i];
1114  (*dQdxMatPtr)[li_mPro[i]][MEquMNodeOffset[i]] += dmQ_dm[i];
1115  (*dQdxMatPtr)[li_hPro[i]][HEquHNodeOffset[i]] += dhQ_dh[i];
1116  }
1117 
1118  return bsuccess;
1119 }
1120 
1121 //-----------------------------------------------------------------------------
1122 // Function : Instance::loadDAEdFdx ()
1123 //
1124 // Purpose : Loads the F-vector contributions for a single
1125 // Neuron 3 instance.
1126 //
1127 // Special Notes : This is an algebraic constaint.
1128 //
1129 // Scope : public
1130 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1131 // Creation Date : 06/10/09
1132 //-----------------------------------------------------------------------------
1134 {
1135  bool bsuccess = true;
1136 
1137  N_LAS_Matrix * dFdxMatPtr = extData.dFdxMatrixPtr;
1138 
1139  (*dFdxMatPtr)[li_Pos][APosEquPosNodeOffset] += dkcl1F_dVin;
1140  (*dFdxMatPtr)[li_Pos][APosEquNextNodeOffset] += dkcl1F_dVs0;
1141 
1142  (*dFdxMatPtr)[li_Neg][ANegEquNegNodeOffset] += dkcl2F_dVout;
1143  (*dFdxMatPtr)[li_Neg][ANegEquLastNodeOffset] += dkcl2F_dVsn;
1144 
1145  for( int i=0; i<nSeg ; i++)
1146  {
1147  (*dFdxMatPtr)[li_Vol[i]][SegVEqnVpreOffset[i]] += segF_dVp[i];
1148  (*dFdxMatPtr)[li_Vol[i]][SegVEqnVsegOffset[i]] += segF_dV[i];
1149  (*dFdxMatPtr)[li_Vol[i]][SegVEqnVnexOffset[i]] += segF_dVn[i];
1150  (*dFdxMatPtr)[li_Vol[i]][SegVEqnNOffset[i]] += segF_dn[i];
1151  (*dFdxMatPtr)[li_Vol[i]][SegVEqnMOffset[i]] += segF_dm[i];
1152  (*dFdxMatPtr)[li_Vol[i]][SegVEqnHOffset[i]] += segF_dh[i];
1153 
1154  (*dFdxMatPtr)[li_nPro[i]][NEquVNodeOffset[i]] += dnF_dV[i];
1155  (*dFdxMatPtr)[li_nPro[i]][NEquNNodeOffset[i]] += dnF_dn[i];
1156  (*dFdxMatPtr)[li_mPro[i]][MEquVNodeOffset[i]] += dmF_dV[i];
1157  (*dFdxMatPtr)[li_mPro[i]][MEquMNodeOffset[i]] += dmF_dm[i];
1158  (*dFdxMatPtr)[li_hPro[i]][HEquVNodeOffset[i]] += dhF_dV[i];
1159  (*dFdxMatPtr)[li_hPro[i]][HEquHNodeOffset[i]] += dhF_dh[i];
1160  }
1161 
1162  return bsuccess;
1163 }
1164 
1165 //-----------------------------------------------------------------------------
1166 // Function : Instance::setIC
1167 // Purpose :
1168 // Special Notes :
1169 // Scope : public
1170 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1171 // Creation Date : 06/10/09
1172 //-----------------------------------------------------------------------------
1174 {
1175  bool bsuccess = true;
1176 
1177  return bsuccess;
1178 }
1179 
1180 //-----------------------------------------------------------------------------
1181 // Function : Instance::varTypes
1182 // Purpose :
1183 // Special Notes :
1184 // Scope : public
1185 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1186 // Creation Date : 06/10/09
1187 //-----------------------------------------------------------------------------
1188 void Instance::varTypes( std::vector<char> & varTypeVec )
1189 {
1190  //varTypeVec.resize(1);
1191  //varTypeVec[0] = 'I';
1192 }
1193 
1194 
1195 //-----------------------------------------------------------------------------
1196 // Function : Model::processParams
1197 // Purpose :
1198 // Special Notes :
1199 // Scope : public
1200 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1201 // Creation Date : 06/10/09
1202 //-----------------------------------------------------------------------------
1204 {
1205  return true;
1206 }
1207 
1208 //----------------------------------------------------------------------------
1209 // Function : Model::processInstanceParams
1210 // Purpose :
1211 // Special Notes :
1212 // Scope : public
1213 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1214 // Creation Date : 06/10/09
1215 //----------------------------------------------------------------------------
1217 {
1218  std::vector<Instance*>::iterator iter;
1219  std::vector<Instance*>::iterator first = instanceContainer.begin();
1220  std::vector<Instance*>::iterator last = instanceContainer.end();
1221 
1222  for (iter=first; iter!=last; ++iter)
1223  {
1224  (*iter)->processParams();
1225  }
1226 
1227  return true;
1228 }
1229 
1230 //-----------------------------------------------------------------------------
1231 // Function : Model::Model
1232 // Purpose : block constructor
1233 // Special Notes :
1234 // Scope : public
1235 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1236 // Creation Date : 06/10/09
1237 //-----------------------------------------------------------------------------
1239  const Configuration & configuration,
1240  const ModelBlock & MB,
1241  const FactoryBlock & factory_block)
1242  : DeviceModel(MB, configuration.getModelParameters(), factory_block),
1243  rInt(0.0),
1244  radius(0.0),
1245  length(0.0),
1246  rIntGiven(false),
1247  radiusGiven(false),
1248  lengthGiven(false)
1249 {
1250 
1251  // Set params to constant default values:
1252  setDefaultParams ();
1253 
1254  // Set params according to .model line and constant defaults from metadata:
1255  setModParams (MB.params);
1256 
1257  // Set any non-constant parameter defaults:
1258  //if (!given("TNOM"))
1259  // tnom = getDeviceOptions().tnom;
1260 
1261  // Calculate any parameters specified as expressions:
1263 
1264  // calculate dependent (ie computed) params and check for errors:
1265 
1266  processParams ();
1267 }
1268 
1269 //-----------------------------------------------------------------------------
1270 // Function : Model::~Model
1271 // Purpose : destructor
1272 // Special Notes :
1273 // Scope : public
1274 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1275 // Creation Date : 06/10/09
1276 //-----------------------------------------------------------------------------
1278 {
1279  std::vector<Instance*>::iterator iter;
1280  std::vector<Instance*>::iterator first = instanceContainer.begin();
1281  std::vector<Instance*>::iterator last = instanceContainer.end();
1282 
1283  for (iter=first; iter!=last; ++iter)
1284  {
1285  delete (*iter);
1286  }
1287 
1288 }
1289 
1290 // additional Declarations
1291 
1292 //-----------------------------------------------------------------------------
1293 // Function : Model::printOutInstances
1294 // Purpose : debugging tool.
1295 // Special Notes :
1296 // Scope : public
1297 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1298 // Creation Date : 06/10/09
1299 //-----------------------------------------------------------------------------
1300 std::ostream &Model::printOutInstances(std::ostream &os) const
1301 {
1302  std::vector<Instance*>::const_iterator iter;
1303  std::vector<Instance*>::const_iterator first = instanceContainer.begin();
1304  std::vector<Instance*>::const_iterator last = instanceContainer.end();
1305 
1306  int i, isize;
1307  isize = instanceContainer.size();
1308 
1309  os << std::endl;
1310  os << "Number of Neuron instances: " << isize << std::endl;
1311  os << " name=\t\tmodelName\tParameters" << std::endl;
1312  for (i=0, iter=first; iter!=last; ++iter, ++i)
1313  {
1314  os << " " << i << ": " << (*iter)->getName() << "\t";
1315  os << getName();
1316  os << std::endl;
1317  }
1318 
1319  os << std::endl;
1320  return os;
1321 }
1322 
1323 //-----------------------------------------------------------------------------
1324 // Function : Model::forEachInstance
1325 // Purpose :
1326 // Special Notes :
1327 // Scope : public
1328 // Creator : David Baur
1329 // Creation Date : 2/4/2014
1330 //-----------------------------------------------------------------------------
1331 /// Apply a device instance "op" to all instances associated with this
1332 /// model
1333 ///
1334 /// @param[in] op Operator to apply to all instances.
1335 ///
1336 ///
1337 void Model::forEachInstance(DeviceInstanceOp &op) const /* override */
1338 {
1339  for (std::vector<Instance *>::const_iterator it = instanceContainer.begin(); it != instanceContainer.end(); ++it)
1340  op(*it);
1341 }
1342 
1343 
1344 Device *Traits::factory(const Configuration &configuration, const FactoryBlock &factory_block)
1345 {
1346 
1347  return new DeviceMaster<Traits>(configuration, factory_block, factory_block.solverState_, factory_block.deviceOptions_);
1348 }
1349 
1351 {
1353  .registerDevice("neuron", 3)
1354  .registerModelType("neuron", 3);
1355 }
1356 
1357 } // namespace Neuron3
1358 } // namespace Device
1359 } // namespace Xyce