Xyce  6.1
N_DEV_Neuron2.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------------
2 // Copyright Notice
3 //
4 // Copyright 2002 Sandia Corporation. Under the terms
5 // of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S.
6 // Government retains certain rights in this software.
7 //
8 // Xyce(TM) Parallel Electrical Simulator
9 // Copyright (C) 2002-2015 Sandia Corporation
10 //
11 // This program is free software: you can redistribute it and/or modify
12 // it under the terms of the GNU General Public License as published by
13 // the Free Software Foundation, either version 3 of the License, or
14 // (at your option) any later version.
15 //
16 // This program is distributed in the hope that it will be useful,
17 // but WITHOUT ANY WARRANTY; without even the implied warranty of
18 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU General Public License
22 // along with this program. If not, see <http://www.gnu.org/licenses/>.
23 //-----------------------------------------------------------------------------
24 
25 //-------------------------------------------------------------------------
26 // Filename : $RCSfile: N_DEV_Neuron2.C,v $
27 //
28 // Purpose :
29 //
30 // Special Notes :
31 //
32 // Creator : Richard Schiek, Electrical and Microsytem Modeling
33 //
34 // Creation Date : 02/28/00
35 //
36 // Revision Information:
37 // ---------------------
38 //
39 // Revision Number: $Revision: 1.49.2.1 $
40 //
41 // Revision Date : $Date: 2015/04/02 18:29:37 $
42 //
43 // Current Owner : $Author: tvrusso $
44 //-------------------------------------------------------------------------
45 
46 #include <Xyce_config.h>
47 
48 
49 // ---------- Standard Includes ----------
50 
51 
52 // ---------- Xyce Includes ----------
53 #include <N_DEV_DeviceOptions.h>
54 #include <N_DEV_DeviceMaster.h>
55 #include <N_DEV_ExternData.h>
56 #include <N_DEV_MatrixLoadData.h>
57 #include <N_DEV_Neuron2.h>
58 #include <N_DEV_SolverState.h>
59 #include <N_DEV_Message.h>
60 #include <N_ERH_ErrorMgr.h>
61 
62 #include <N_LAS_Vector.h>
63 #include <N_LAS_Matrix.h>
64 #include <N_UTL_FeatureTest.h>
65 
66 namespace Xyce {
67 namespace Device {
68 
69 namespace Neuron2 {
70 
72 {
73 }
74 
76 {
77  p.addPar ("CMEM",0.0,&Neuron2::Model::cMem)
78  .setGivenMember(&Neuron2::Model::cMemGiven)
79  .setUnit(U_FARAD)
80  .setCategory(CAT_NONE)
81  .setDescription("Membrane capacitance");
82 
83  p.addPar ("GMEM",0.0,&Neuron2::Model::gMem)
84  .setGivenMember(&Neuron2::Model::gMemGiven)
85  .setUnit(U_OHMM1)
86  .setCategory(CAT_NONE)
87  .setDescription("Membrane conductance");
88 
89  p.addPar ("VREST",0.0,&Neuron2::Model::vRest)
90  .setGivenMember(&Neuron2::Model::vRestGiven)
91  .setUnit(U_VOLT)
92  .setCategory(CAT_NONE)
93  .setDescription("Resting potential");
94 
95  p.addPar ("EK",0.0,&Neuron2::Model::eK)
96  .setGivenMember(&Neuron2::Model::eKGiven)
97  .setUnit(U_VOLT)
98  .setCategory(CAT_NONE)
99  .setDescription("Potassium resting potential");
100 
101  p.addPar ("GK",0.0,&Neuron2::Model::gK)
102  .setGivenMember(&Neuron2::Model::gKGiven)
103  .setUnit(U_OHMM1)
104  .setCategory(CAT_NONE)
105  .setDescription("Potassium base conductance");
106 
107  p.addPar ("ENA",0.0,&Neuron2::Model::eNa)
108  .setGivenMember(&Neuron2::Model::eNaGiven)
109  .setUnit(U_VOLT)
110  .setCategory(CAT_NONE)
111  .setDescription("Sodium resting potential");
112 
113  p.addPar ("GNA",0.0,&Neuron2::Model::gNa)
114  .setGivenMember(&Neuron2::Model::gNaGiven)
115  .setUnit(U_OHMM1)
116  .setCategory(CAT_NONE)
117  .setDescription("Sodium base conductance");
118 
119  p.addPar ("EA",0.0,&Neuron2::Model::eA)
120  .setGivenMember(&Neuron2::Model::eAGiven)
121  .setUnit(U_CM)
122  .setCategory(CAT_NONE)
123  .setDescription("a-current rest potential");
124 
125  p.addPar ("GA",0.0,&Neuron2::Model::gA)
126  .setGivenMember(&Neuron2::Model::gAGiven)
127  .setUnit(U_NONE)
128  .setCategory(CAT_NONE)
129  .setDescription("a-current base conductance");
130 
131  p.addPar ("ECA",0.0,&Neuron2::Model::eCa)
132  .setGivenMember(&Neuron2::Model::eCaGiven)
133  .setUnit(U_NONE)
134  .setCategory(CAT_NONE)
135  .setDescription("Calcium rest potential");
136 
137  p.addPar ("GCA",0.0,&Neuron2::Model::gCa)
138  .setGivenMember(&Neuron2::Model::gCaGiven)
139  .setUnit(U_OHMM1)
140  .setCategory(CAT_NONE)
141  .setDescription("Calcium base conductance");
142 
143  p.addPar ("EKCA",0.0,&Neuron2::Model::eKCa)
144  .setGivenMember(&Neuron2::Model::eKCaGiven)
145  .setUnit(U_OHMM1)
146  .setCategory(CAT_NONE)
147  .setDescription("Potassium-calcium rest potential");
148 
149  p.addPar ("GKCA",0.0,&Neuron2::Model::gKCa)
150  .setGivenMember(&Neuron2::Model::gKCaGiven)
151  .setUnit(U_OHMM1)
152  .setCategory(CAT_NONE)
153  .setDescription("Potassium-calcium base conductance");
154 
155  p.addPar ("CAINIT",0.0,&Neuron2::Model::CaInit)
156  .setGivenMember(&Neuron2::Model::CaInitGiven)
157  .setUnit(U_OHMM1)
158  .setCategory(CAT_NONE)
159  .setDescription("initial intra-cellular calcium concentration");
160 
161  p.addPar ("CAGAMMA",0.0,&Neuron2::Model::CaGamma)
162  .setGivenMember(&Neuron2::Model::CaGammaGiven)
163  .setUnit(U_OHMM1)
164  .setCategory(CAT_NONE)
165  .setDescription("calcium current to concentration multiplier");
166 
167  p.addPar ("CATAU",0.0,&Neuron2::Model::CaTau)
168  .setGivenMember(&Neuron2::Model::CaTauGiven)
169  .setUnit(U_OHMM1)
170  .setCategory(CAT_NONE)
171  .setDescription("calcium removal time constant");
172 }
173 
174 //
175 // static class member inits
176 //
177 std::vector< std::vector<int> > Instance::jacStamp;
178 
179 // Class Instance
180 //-----------------------------------------------------------------------------
181 // Function : Instance::processParams
182 // Purpose :
183 // Special Notes :
184 // Scope : public
185 // Creator : Richard Schiek, Electrical and Microsytem Modeling
186 // Creation Date : 01/02/08
187 //-----------------------------------------------------------------------------
189 {
190  // If there are any time dependent parameters, set their values at for
191  // the current time.
192 
193  // now set the temperature related stuff.
194  //updateTemperature(temp);
195 
196  return true;
197 }
198 
199 //-----------------------------------------------------------------------------
200 // Function : Instance::updateTemperature
201 // Purpose :
202 // Special Notes :
203 // Scope : public
204 // Creator : Richard Schiek, Electrical and Microsytem Modeling
205 // Creation Date : 01/02/08
206 //-----------------------------------------------------------------------------
207 bool Instance::updateTemperature ( const double & temp)
208 {
209  bool bsuccess = true;
210  return bsuccess;
211 }
212 
213 //-----------------------------------------------------------------------------
214 // Function : Instance::Instance
215 // Purpose : constructor
216 // Special Notes :
217 // Scope : public
218 // Creator : Richard Schiek, Electrical and Microsytem Modeling
219 // Creation Date : 01/02/08
220 //-----------------------------------------------------------------------------
222  const Configuration & configuration,
223  const InstanceBlock & IB,
224  Model & Miter,
225  const FactoryBlock & factory_block)
226  : DeviceInstance(IB, configuration.getInstanceParameters(), factory_block),
227  model_(Miter),
228  kcl1Fvalue(0.0),
229  kcl1Qvalue(0.0),
230  kcl2Fvalue(0.0),
231  kcl2Qvalue(0.0),
232  nEquFvalue(0.0),
233  nEquQvalue(0.0),
234  mEquFvalue(0.0),
235  mEquQvalue(0.0),
236  hEquFvalue(0.0),
237  hEquQvalue(0.0),
238  dkcl1F_dV1(0.0),
239  dkcl1F_dV2(0.0),
240  dkcl1F_dn(0.0),
241  dkcl1F_dm(0.0),
242  dkcl1F_dh(0.0),
243  dkcl1Q_dV1(0.0),
244  dkcl1Q_dV2(0.0),
245  dkcl2F_dV1(0.0),
246  dkcl2F_dV2(0.0),
247  dkcl2F_dn(0.0),
248  dkcl2F_dm(0.0),
249  dkcl2F_dh(0.0),
250  dkcl2Q_dV1(0.0),
251  dkcl2Q_dV2(0.0),
252  dnF_dV1(0.0),
253  dnF_dn(0.0),
254  dnQ_dn(0.0),
255  dmF_dV1(0.0),
256  dmF_dm(0.0),
257  dmQ_dm(0.0),
258  dhF_dV1(0.0),
259  dhF_dh(0.0),
260  dhQ_dh(0.0)
261 {
262  numExtVars = 2;
263  numIntVars = 9;
264  numStateVars = 2;
265 
266  // set up jacStamp
267  if( jacStamp.empty() )
268  {
269  jacStamp.resize(11);
270  jacStamp[0].resize(10); // kcl 1
271  jacStamp[0][0] = 0;
272  jacStamp[0][1] = 1;
273  jacStamp[0][2] = 2;
274  jacStamp[0][3] = 3;
275  jacStamp[0][4] = 4;
276  jacStamp[0][5] = 5;
277  jacStamp[0][6] = 6;
278  jacStamp[0][7] = 7;
279  jacStamp[0][8] = 8;
280  jacStamp[0][9] = 9;
281  jacStamp[1].resize(10); // kcl 2
282  jacStamp[1][0] = 0;
283  jacStamp[1][1] = 1;
284  jacStamp[1][2] = 2;
285  jacStamp[1][3] = 3;
286  jacStamp[1][4] = 4;
287  jacStamp[1][5] = 5;
288  jacStamp[1][6] = 6;
289  jacStamp[1][7] = 7;
290  jacStamp[1][8] = 8;
291  jacStamp[1][9] = 9;
292  jacStamp[2].resize(2); // n equation
293  jacStamp[2][0] = 0;
294  jacStamp[2][1] = 2;
295  jacStamp[3].resize(2); // m equation
296  jacStamp[3][0] = 0;
297  jacStamp[3][1] = 3;
298  jacStamp[4].resize(2); // h equation
299  jacStamp[4][0] = 0;
300  jacStamp[4][1] = 4;
301  jacStamp[5].resize(2); // a equation
302  jacStamp[5][0] = 0;
303  jacStamp[5][1] = 5;
304  jacStamp[6].resize(2); // b equation
305  jacStamp[6][0] = 0;
306  jacStamp[6][1] = 6;
307  jacStamp[7].resize(2); // M equation
308  jacStamp[7][0] = 0;
309  jacStamp[7][1] = 7;
310  jacStamp[8].resize(2); // H equation
311  jacStamp[8][0] = 0;
312  jacStamp[8][1] = 8;
313  jacStamp[9].resize(3); // c equation
314  jacStamp[9][0] = 0;
315  jacStamp[9][1] = 9;
316  jacStamp[9][2] = 10;
317  jacStamp[10].resize(5); // Ca equation
318  jacStamp[10][0] = 0;
319  jacStamp[10][1] = 1;
320  jacStamp[10][2] = 7;
321  jacStamp[10][3] = 8;
322  jacStamp[10][4] = 10;
323  }
324 
325  // Set params to constant default values:
326  setDefaultParams ();
327 
328  // Set params according to instance line and constant defaults from metadata:
329  // there are no params for this device as yet, so calling setParams causes
330  // the code to exit.
331  // setParams (IB.params);
332 
333  // Set any non-constant parameter defaults:
334 
335  //if (!given("TEMP"))
336  // temp = getDeviceOptions().temp.dVal();
337 
338  // Calculate any parameters specified as expressions:
340 
341  // calculate dependent (ie computed) params and check for errors:
342  processParams ();
343 
344 }
345 
346 //-----------------------------------------------------------------------------
347 // Function : Instance::~Instance
348 // Purpose : destructor
349 // Special Notes :
350 // Scope : public
351 // Creator : Richard Schiek, Electrical and Microsytem Modeling
352 // Creation Date : 01/02/08
353 //-----------------------------------------------------------------------------
355 {
356 }
357 
358 //-----------------------------------------------------------------------------
359 // Function : Instance::registerLIDs
360 // Purpose :
361 // Special Notes :
362 // Scope : public
363 // Creator : Richard Schiek, Electrical and Microsytem Modeling
364 // Creation Date : 01/02/08
365 //-----------------------------------------------------------------------------
366 void Instance::registerLIDs(const std::vector<int> & intLIDVecRef,
367  const std::vector<int> & extLIDVecRef)
368 {
369  AssertLIDs(intLIDVecRef.size() == numIntVars);
370  AssertLIDs(extLIDVecRef.size() == numExtVars);
371 
372  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
373  {
374  Xyce::dout() << std::endl << section_divider << std::endl;
375  Xyce::dout() << " Instance::registerLIDs" << std::endl;
376  Xyce::dout() << " name = " << getName() << std::endl;
377  }
378 
379  // copy over the global ID lists.
380  intLIDVec = intLIDVecRef;
381  extLIDVec = extLIDVecRef;
382 
383  li_Pos = extLIDVec[0];
384  li_Neg = extLIDVec[1];
385  li_nPro = intLIDVec[0];
386  li_mPro = intLIDVec[1];
387  li_hPro = intLIDVec[2];
388  li_aPro = intLIDVec[3];
389  li_bPro = intLIDVec[4];
390  li_M_Pro = intLIDVec[5];
391  li_H_Pro = intLIDVec[6];
392  li_cPro = intLIDVec[7];
393  li_CaPro = intLIDVec[8];
394 
395  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) )
396  {
397  Xyce::dout() << " li_Pos = " << li_Pos << std::endl
398  << " li_Neg = " << li_Neg << std::endl
399  << " li_nPro = " << li_nPro << std::endl
400  << " li_mPro = " << li_mPro << std::endl
401  << " li_hPro = " << li_hPro << std::endl
402  << " li_aPro = " << li_aPro << std::endl
403  << " li_bPro = " << li_bPro << std::endl
404  << " li_M_Pro = " << li_M_Pro << std::endl
405  << " li_H_Pro = " << li_H_Pro << std::endl
406  << " li_cPro = " << li_cPro << std::endl
407  << " li_CaPro = " << li_CaPro << std::endl;
408  }
409 
410  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) )
411  {
412  Xyce::dout() << section_divider << std::endl;
413  }
414 }
415 
416 //-----------------------------------------------------------------------------
417 // Function : Instance::loadNodeSymbols
418 // Purpose :
419 // Special Notes :
420 // Scope : public
421 // Creator : Eric R. Keiter, SNL, Parallel Computational Sciences
422 // Creation Date : 05/13/05
423 //-----------------------------------------------------------------------------
424 void Instance::loadNodeSymbols(Util::SymbolTable &symbol_table) const
425 {
426  addInternalNode(symbol_table, li_nPro, getName(), "N");
427  addInternalNode(symbol_table, li_mPro, getName(), "M");
428  addInternalNode(symbol_table, li_hPro, getName(), "H");
429  addInternalNode(symbol_table, li_aPro, getName(), "A");
430  addInternalNode(symbol_table, li_bPro, getName(), "B");
431  addInternalNode(symbol_table, li_M_Pro, getName(), "M_");
432  addInternalNode(symbol_table, li_H_Pro, getName(), "H_");
433  addInternalNode(symbol_table, li_cPro, getName(), "C");
434  addInternalNode(symbol_table, li_CaPro, getName(), "Ca");
435 
436 }
437 
438 //-----------------------------------------------------------------------------
439 // Function : Instance::registerStateLIDs
440 // Purpose :
441 // Special Notes :
442 // Scope : public
443 // Creator : Richard Schiek, Electrical and Microsytem Modeling
444 // Creation Date : 01/02/08
445 //-----------------------------------------------------------------------------
446 void Instance::registerStateLIDs( const std::vector<int> & staLIDVecRef )
447 {
448  AssertLIDs(staLIDVecRef.size() == numStateVars);
449 
450  // copy over the global ID lists.
451  staLIDVec = staLIDVecRef;
452 
453  li_KCurrentState = staLIDVec[0];
454  li_NaCurrentState = staLIDVec[1];
455 
456 }
457 
458 //-----------------------------------------------------------------------------
459 // Function : Instance::jacobianStamp
460 // Purpose :
461 // Special Notes :
462 // Scope : public
463 // Creator : Richard Schiek, Electrical and Microsytem Modeling
464 // Creation Date : 01/02/08
465 //-----------------------------------------------------------------------------
466 const std::vector< std::vector<int> > & Instance::jacobianStamp() const
467 {
468  return jacStamp;
469 }
470 
471 //-----------------------------------------------------------------------------
472 // Function : Instance::registerJacLIDs
473 // Purpose :
474 // Special Notes :
475 // Scope : public
476 // Creator : Richard Schiek, Electrical and Microsytem Modeling
477 // Creation Date : 01/02/08
478 //-----------------------------------------------------------------------------
479 void Instance::registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec )
480 {
481  DeviceInstance::registerJacLIDs( jacLIDVec );
482 
483  APosEquPosNodeOffset = jacLIDVec[0][0];
484  APosEquNegNodeOffset = jacLIDVec[0][1];
485  APosEquNNodeOffset = jacLIDVec[0][2];
486  APosEquMNodeOffset = jacLIDVec[0][3];
487  APosEquHNodeOffset = jacLIDVec[0][4];
488  APosEquANodeOffset = jacLIDVec[0][5];
489  APosEquBNodeOffset = jacLIDVec[0][6];
490  APosEquM_NodeOffset = jacLIDVec[0][7];
491  APosEquH_NodeOffset = jacLIDVec[0][8];
492  APosEquCNodeOffset = jacLIDVec[0][9];
493 
494  ANegEquPosNodeOffset = jacLIDVec[1][0];
495  ANegEquNegNodeOffset = jacLIDVec[1][1];
496  ANegEquNNodeOffset = jacLIDVec[1][2];
497  ANegEquMNodeOffset = jacLIDVec[1][3];
498  ANegEquHNodeOffset = jacLIDVec[1][4];
499  ANegEquANodeOffset = jacLIDVec[1][5];
500  ANegEquBNodeOffset = jacLIDVec[1][6];
501  ANegEquM_NodeOffset = jacLIDVec[1][7];
502  ANegEquH_NodeOffset = jacLIDVec[1][8];
503  ANegEquCNodeOffset = jacLIDVec[1][9];
504 
505  ANEquPosNodeOffset = jacLIDVec[2][0];
506  ANEquNNodeOffset = jacLIDVec[2][1];
507 
508  AMEquPosNodeOffset = jacLIDVec[3][0];
509  AMEquMNodeOffset = jacLIDVec[3][1];
510 
511  AHEquPosNodeOffset = jacLIDVec[4][0];
512  AHEquHNodeOffset = jacLIDVec[4][1];
513 
514  AAEquPosNodeOffset = jacLIDVec[5][0];
515  AAEquANodeOffset = jacLIDVec[5][1];
516 
517  ABEquPosNodeOffset = jacLIDVec[6][0];
518  ABEquBNodeOffset = jacLIDVec[6][1];
519 
520  AM_EquPosNodeOffset = jacLIDVec[7][0];
521  AM_EquM_NodeOffset = jacLIDVec[7][1];
522 
523  AH_EquPosNodeOffset = jacLIDVec[8][0];
524  AH_EquH_NodeOffset = jacLIDVec[8][1];
525 
526  ACEquPosNodeOffset = jacLIDVec[9][0];
527  ACEquCNodeOffset = jacLIDVec[9][1];
528  ACEquCaNodeOffset = jacLIDVec[9][2];
529 
530  ACaEquPosNodeOffset = jacLIDVec[10][0];
531  ACaEquNegNodeOffset = jacLIDVec[10][1];
532  ACaEquM_NodeOffset = jacLIDVec[10][2];
533  ACaEquH_NodeOffset = jacLIDVec[10][3];
534  ACaEquCaNodeOffset = jacLIDVec[10][4];
535 
536 }
537 
538 
539 //-----------------------------------------------------------------------------
540 // Function : Instance::updateIntermediateVars
541 // Purpose :
542 // Special Notes :
543 // Scope : public
544 // Creator : Richard Schiek, Electrical and Microsytem Modeling
545 // Creation Date : 01/02/08
546 //-----------------------------------------------------------------------------
548 {
549  //Xyce::dout() << "Instance::updateIntermediateVars()" << std::endl;
550 
551  bool bsuccess = true;
552 
553  // here we take the current solutions for V1, V2, n, m and h
554  // and use those to calculate all the terms needed for the next
555  // load cycle (F, Q, dFdX, dQdX)
556  Linear::Vector * solVectorPtr = extData.nextSolVectorPtr;
557 
558  // use suffix "now" to to clarify that this for the latest solution
559  double v1Now = (*solVectorPtr)[li_Pos];
560  double v2Now = (*solVectorPtr)[li_Neg];
561  double nNow = (*solVectorPtr)[li_nPro];
562  double mNow = (*solVectorPtr)[li_mPro];
563  double hNow = (*solVectorPtr)[li_hPro];
564  double aNow = (*solVectorPtr)[li_aPro];
565  double bNow = (*solVectorPtr)[li_bPro];
566  double M_Now = (*solVectorPtr)[li_M_Pro];
567  double H_Now = (*solVectorPtr)[li_H_Pro];
568  double cNow = (*solVectorPtr)[li_cPro];
569  double CaNow = (*solVectorPtr)[li_CaPro];
570 
571  // get function and derivative values
572  // independent variables
573  // use scoping to avoid lots of similar variable names
574  {
575  Sacado::Fad::SFad<double,10> v1Var( 10, 0, v1Now );
576  Sacado::Fad::SFad<double,10> v2Var( 10, 1, v2Now );
577  Sacado::Fad::SFad<double,10> nVar( 10, 2, nNow );
578  Sacado::Fad::SFad<double,10> mVar( 10, 3, mNow );
579  Sacado::Fad::SFad<double,10> hVar( 10, 4, hNow );
580  Sacado::Fad::SFad<double,10> aVar( 10, 5, aNow );
581  Sacado::Fad::SFad<double,10> bVar( 10, 6, bNow );
582  Sacado::Fad::SFad<double,10> M_Var( 10, 7, M_Now );
583  Sacado::Fad::SFad<double,10> H_Var( 10, 8, H_Now );
584  Sacado::Fad::SFad<double,10> cVar( 10, 9, cNow );
585 
586  // parameters from the model that we'll need.
587  Sacado::Fad::SFad<double,10> gMemVar( model_.gMem );
588  Sacado::Fad::SFad<double,10> vRestVar( model_.vRest );
589  Sacado::Fad::SFad<double,10> gKVar( model_.gK );
590  Sacado::Fad::SFad<double,10> eKVar( model_.eK );
591  Sacado::Fad::SFad<double,10> gNaVar( model_.gNa );
592  Sacado::Fad::SFad<double,10> eNaVar( model_.eNa );
593  Sacado::Fad::SFad<double,10> gAVar( model_.gA );
594  Sacado::Fad::SFad<double,10> eAVar( model_.eA );
595  Sacado::Fad::SFad<double,10> gCaVar( model_.gCa );
596  Sacado::Fad::SFad<double,10> eCaVar( model_.eCa );
597  Sacado::Fad::SFad<double,10> gKCaVar( model_.gKCa );
598  Sacado::Fad::SFad<double,10> CaInitVar( model_.CaInit );
599  Sacado::Fad::SFad<double,10> CaGammaVar( model_.CaGamma );
600  Sacado::Fad::SFad<double,10> CaTauVar( model_.CaTau );
601 
602  // compute the vaud and derivative terms for KCL 1 F
603  Sacado::Fad::SFad<double,10> resultFad;
604  resultFad = kcl1EquF( v1Var, v2Var, nVar, mVar, hVar, aVar, bVar, M_Var, H_Var, cVar, gMemVar, vRestVar, gKVar, eKVar, gNaVar, eNaVar, gAVar, eAVar, gCaVar, eCaVar, gKCaVar);
605  kcl1Fvalue = resultFad.val();
606  dkcl1F_dV1 = resultFad.dx(0);
607  dkcl1F_dV2 = resultFad.dx(1);
608  dkcl1F_dn = resultFad.dx(2);
609  dkcl1F_dm = resultFad.dx(3);
610  dkcl1F_dh = resultFad.dx(4);
611  dkcl1F_da = resultFad.dx(5);
612  dkcl1F_db = resultFad.dx(6);
613  dkcl1F_dM = resultFad.dx(7);
614  dkcl1F_dH = resultFad.dx(8);
615  dkcl1F_dc = resultFad.dx(9);
616 
617  // compute the vaud and derivative terms for KCL 2 F
618  resultFad = kcl2EquF( v1Var, v2Var, nVar, mVar, hVar, aVar, bVar, M_Var, H_Var, cVar, gMemVar, vRestVar, gKVar, eKVar, gNaVar, eNaVar, gAVar, eAVar, gCaVar, eCaVar, gKCaVar);
619  kcl2Fvalue = resultFad.val();
620  dkcl2F_dV1 = resultFad.dx(0);
621  dkcl2F_dV2 = resultFad.dx(1);
622  dkcl2F_dn = resultFad.dx(2);
623  dkcl2F_dm = resultFad.dx(3);
624  dkcl2F_dh = resultFad.dx(4);
625  dkcl2F_da = resultFad.dx(5);
626  dkcl2F_db = resultFad.dx(6);
627  dkcl2F_dM = resultFad.dx(7);
628  dkcl2F_dH = resultFad.dx(8);
629  dkcl2F_dc = resultFad.dx(9);
630  }
631 
632  {
633  Sacado::Fad::SFad<double,2> v1Var( 2, 0, v1Now );
634  Sacado::Fad::SFad<double,2> v2Var( 2, 1, v2Now );
635 
636  // parameters
637  Sacado::Fad::SFad<double,2> cMemVar( model_.cMem );
638 
639  Sacado::Fad::SFad<double,2> resultFad;
640  resultFad = kcl1EquQ( v1Var, v2Var, cMemVar );
641  kcl1Qvalue = resultFad.val();
642  dkcl1Q_dV1 = resultFad.dx(0);
643  dkcl1Q_dV2 = resultFad.dx(1);
644 
645  resultFad = kcl2EquQ( v1Var, v2Var, cMemVar );
646  kcl2Qvalue = resultFad.val();
647  dkcl2Q_dV1 = resultFad.dx(0);
648  dkcl2Q_dV2 = resultFad.dx(1);
649  }
650 
651  // n - equation
652  {
653  Sacado::Fad::SFad<double,2> v1Var( 2, 0, v1Now );
654  Sacado::Fad::SFad<double,2> nVar( 2, 1, nNow );
655  // parameter
656  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
657 
658  Sacado::Fad::SFad<double,2> resultFad = nEquF( v1Var, nVar, vRestVar);
659  nEquFvalue = resultFad.val();
660  dnF_dV1 = resultFad.dx(0);
661  dnF_dn = resultFad.dx(1);
662  }
663 
664  {
665  Sacado::Fad::SFad<double,1> nVar( 1, 0, nNow );
666 
667  Sacado::Fad::SFad<double,1> resultFad = nEquQ( nVar );
668  nEquQvalue = resultFad.val();
669  dnQ_dn = resultFad.dx(0);
670  }
671 
672  // m - equation
673  {
674  Sacado::Fad::SFad<double,2> v1Var( 2, 0, v1Now );
675  Sacado::Fad::SFad<double,2> mVar( 2, 1, mNow );
676  // parameter
677  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
678 
679  Sacado::Fad::SFad<double,2> resultFad = mEquF( v1Var, mVar, vRestVar );
680  mEquFvalue = resultFad.val();
681  dmF_dV1 = resultFad.dx(0);
682  dmF_dm = resultFad.dx(1);
683  }
684  {
685  Sacado::Fad::SFad<double,1> mVar( 1, 0, mNow );
686 
687  Sacado::Fad::SFad<double,1> resultFad = mEquQ( mVar );
688  mEquQvalue = resultFad.val();
689  dmQ_dm = resultFad.dx(0);
690  }
691 
692  // h - equation
693  {
694  Sacado::Fad::SFad<double,2> v1Var( 2, 0, v1Now );
695  Sacado::Fad::SFad<double,2> hVar( 2, 1, hNow );
696  // parameter
697  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
698 
699  Sacado::Fad::SFad<double,2> resultFad = hEquF( v1Var, hVar, vRestVar );
700  hEquFvalue = resultFad.val();
701  dhF_dV1 = resultFad.dx(0);
702  dhF_dh = resultFad.dx(1);
703  }
704  {
705  Sacado::Fad::SFad<double,1> hVar( 1, 0, hNow );
706 
707  Sacado::Fad::SFad<double,1> resultFad = hEquQ( hVar );
708  hEquQvalue = resultFad.val();
709  dhQ_dh = resultFad.dx(0);
710  }
711 
712  // a - equation
713  {
714  Sacado::Fad::SFad<double,2> v1Var( 2, 0, v1Now );
715  Sacado::Fad::SFad<double,2> aVar( 2, 1, aNow );
716  // parameter
717  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
718 
719  Sacado::Fad::SFad<double,2> resultFad = aEquF( v1Var, aVar, vRestVar );
720  aEquFvalue = resultFad.val();
721  daF_dV1 = resultFad.dx(0);
722  daF_da = resultFad.dx(1);
723  }
724  {
725  Sacado::Fad::SFad<double,1> aVar( 1, 0, aNow );
726 
727  Sacado::Fad::SFad<double,1> resultFad = aEquQ( aVar );
728  aEquQvalue = resultFad.val();
729  daQ_da = resultFad.dx(0);
730  }
731 
732  // b - equation
733  {
734  Sacado::Fad::SFad<double,2> v1Var( 2, 0, v1Now );
735  Sacado::Fad::SFad<double,2> bVar( 2, 1, bNow );
736  // parameter
737  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
738 
739  Sacado::Fad::SFad<double,2> resultFad = bEquF( v1Var, bVar, vRestVar );
740  bEquFvalue = resultFad.val();
741  dbF_dV1 = resultFad.dx(0);
742  dbF_db = resultFad.dx(1);
743  }
744  {
745  Sacado::Fad::SFad<double,1> bVar( 1, 0, bNow );
746 
747  Sacado::Fad::SFad<double,1> resultFad = aEquQ( bVar );
748  bEquQvalue = resultFad.val();
749  dbQ_db = resultFad.dx(0);
750  }
751 
752  // M - equation
753  {
754  Sacado::Fad::SFad<double,2> v1Var( 2, 0, v1Now );
755  Sacado::Fad::SFad<double,2> M_Var( 2, 1, M_Now );
756  // parameter
757  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
758 
759  Sacado::Fad::SFad<double,2> resultFad = M_EquF( v1Var, M_Var, vRestVar );
760  M_EquFvalue = resultFad.val();
761  dMF_dV1 = resultFad.dx(0);
762  dMF_dM = resultFad.dx(1);
763  }
764  {
765  Sacado::Fad::SFad<double,1> M_Var( 1, 0, M_Now );
766 
767  Sacado::Fad::SFad<double,1> resultFad = aEquQ( M_Var );
768  M_EquQvalue = resultFad.val();
769  dMQ_dM = resultFad.dx(0);
770  }
771 
772  // H - equation
773  {
774  Sacado::Fad::SFad<double,2> v1Var( 2, 0, v1Now );
775  Sacado::Fad::SFad<double,2> H_Var( 2, 1, H_Now );
776  // parameter
777  Sacado::Fad::SFad<double,2> vRestVar( model_.vRest );
778 
779  Sacado::Fad::SFad<double,2> resultFad = H_EquF( v1Var, H_Var, vRestVar );
780  H_EquFvalue = resultFad.val();
781  dHF_dV1 = resultFad.dx(0);
782  dHF_dH = resultFad.dx(1);
783  }
784  {
785  Sacado::Fad::SFad<double,1> H_Var( 1, 0, H_Now );
786 
787  Sacado::Fad::SFad<double,1> resultFad = H_EquQ( H_Var );
788  H_EquQvalue = resultFad.val();
789  dHQ_dH = resultFad.dx(0);
790  }
791 
792  // c - equation
793  {
794  Sacado::Fad::SFad<double,3> v1Var( 3, 0, v1Now );
795  Sacado::Fad::SFad<double,3> cVar( 3, 1, cNow );
796  Sacado::Fad::SFad<double,3> CaVar( 3, 2, CaNow );
797  // parameter
798  Sacado::Fad::SFad<double,3> vRestVar( model_.vRest );
799 
800  Sacado::Fad::SFad<double,3> resultFad = C_EquF( v1Var, cVar, CaVar, vRestVar );
801  cEquFvalue = resultFad.val();
802  dcF_dV1 = resultFad.dx(0);
803  dcF_dc = resultFad.dx(1);
804  dcF_dCa = resultFad.dx(2);
805  }
806  {
807  Sacado::Fad::SFad<double,1> cVar( 1, 0, cNow );
808 
809  Sacado::Fad::SFad<double,1> resultFad = C_EquQ( cVar );
810  cEquQvalue = resultFad.val();
811  dcQ_dc = resultFad.dx(0);
812  }
813 
814  // Ca - equation
815  {
816  Sacado::Fad::SFad<double,5> v1Var( 5, 0, v1Now );
817  Sacado::Fad::SFad<double,5> v2Var( 5, 1, v2Now );
818  Sacado::Fad::SFad<double,5> M_Var( 5, 2, M_Now );
819  Sacado::Fad::SFad<double,5> H_Var( 5, 3, H_Now );
820  Sacado::Fad::SFad<double,5> CaVar( 5, 4, CaNow );
821 
822  // parameter
823  Sacado::Fad::SFad<double,5> gCaVar( model_.gCa );
824  Sacado::Fad::SFad<double,5> eCaVar( model_.gCa );
825  Sacado::Fad::SFad<double,5> CaGammaVar( model_.CaGamma );
826  Sacado::Fad::SFad<double,5> CaTauVar( model_.CaTau );
827 
828  Sacado::Fad::SFad<double,5> resultFad = Ca_EquF( v1Var, v2Var, M_Var, H_Var, CaVar, gCaVar, eCaVar, CaGammaVar, CaTauVar );
829  CaEquFvalue = resultFad.val();
830  dCaF_dV1 = resultFad.dx(0);
831  dCaF_dV2 = resultFad.dx(1);
832  dCaF_dM = resultFad.dx(2);
833  dCaF_dH = resultFad.dx(3);
834  dCaF_dCa = resultFad.dx(4);
835  }
836  {
837  Sacado::Fad::SFad<double,1> CaVar( 1, 0, CaNow );
838 
839  Sacado::Fad::SFad<double,1> resultFad = Ca_EquQ( CaVar );
840  CaEquQvalue = resultFad.val();
841  dCaQ_dCa = resultFad.dx(0);
842  }
843  // calculate potassium current
844  potassiumCurrent = model_.gK * pow(nNow, 4.0) * (v1Now - v2Now - model_.eK);
845 
846  // calculate sodium current
847  sodiumCurrent = model_.gNa * pow(mNow, 3.0) * hNow * (v1Now - v2Now - model_.eNa);
848 
849 
850  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
851  {
852  Xyce::dout() << "Instance::updateIntermediateVars()" << std::endl
853  << "v1 = " << v1Now << std::endl
854  << "v2 = " << v2Now << std::endl
855  << "nNow = " << nNow << std::endl
856  << "mNow = " << mNow << std::endl
857  << "hNow = " << hNow << std::endl
858  << "aNow = " << aNow << std::endl
859  << "bNow = " << bNow << std::endl
860  << "M_Now = " << M_Now << std::endl
861  << "H_Now = " << H_Now << std::endl
862  << "cNow = " << cNow << std::endl
863  << "CaNow = " << CaNow << std::endl
864  << "kcl1Fvalue = " << kcl1Fvalue << std::endl
865  << "dkcl1F_dV1 = " << dkcl1F_dV1 << std::endl
866  << "dkcl1F_dV2 = " << dkcl1F_dV2 << std::endl
867  << "dkcl1F_dn = " << dkcl1F_dn << std::endl
868  << "dkcl1F_dm = " << dkcl1F_dm << std::endl
869  << "dkcl1F_dh = " << dkcl1F_dh << std::endl
870  << "kcl2Fvalue = " << kcl2Fvalue << std::endl
871  << "dkcl2F_dV1 = " << dkcl2F_dV1 << std::endl
872  << "dkcl2F_dV2 = " << dkcl2F_dV2 << std::endl
873  << "dkcl2F_dn = " << dkcl2F_dn << std::endl
874  << "dkcl2F_dm = " << dkcl2F_dm << std::endl
875  << "dkcl2F_dh = " << dkcl2F_dh << std::endl
876  << "alphaN = " << alphaN<double>( v1Now ) << std::endl
877  << "betaN = " << betaN<double>( v1Now ) << std::endl
878  << "nEquFvalue = " << nEquFvalue << std::endl
879  << "dnF_dV1 = " << dnF_dV1 << std::endl
880  << "dnF_dn = " << dnF_dn << std::endl
881  << "nEquQvalue = " << nEquQvalue << std::endl
882  << "dnQ_dn = " << dnQ_dn << std::endl
883  << "alphaM = " << alphaM<double>( v1Now ) << std::endl
884  << "betaM = " << betaM<double>( v1Now ) << std::endl
885  << "mEquFvalue = " << mEquFvalue << std::endl
886  << "dmF_dV1 = " << dmF_dV1 << std::endl
887  << "dmF_dm = " << dmF_dm << std::endl
888  << "mEquQvalue = " << mEquQvalue << std::endl
889  << "dmQ_dm = " << dmQ_dm << std::endl
890  << "alphaH = " << alphaH<double>( v1Now ) << std::endl
891  << "betaH = " << betaH<double>( v1Now ) << std::endl
892  << "hEquFvalue = " << hEquFvalue << std::endl
893  << "dhF_dV1 = " << dhF_dV1 << std::endl
894  << "dhF_dh = " << dhF_dh << std::endl
895  << "hEquQvalue = " << hEquQvalue << std::endl
896  << "dhQ_dh = " << dhQ_dh << std::endl
897 
898  << "aInf = " << aInf<double>( v1Now ) << std::endl
899  << "aTau = " << aTau<double>( v1Now ) << std::endl
900  << "aEquFvalue = " << aEquFvalue << std::endl
901  << "daF_dV1 = " << daF_dV1 << std::endl
902  << "daF_da = " << daF_da << std::endl
903  << "aEquQvalue = " << aEquQvalue << std::endl
904  << "daQ_da = " << daQ_da << std::endl
905 
906  << "bInf = " << bInf<double>( v1Now ) << std::endl
907  << "bTau = " << bTau<double>( v1Now ) << std::endl
908  << "bEquFvalue = " << bEquFvalue << std::endl
909  << "dbF_dV1 = " << dbF_dV1 << std::endl
910  << "dbF_db = " << dbF_db << std::endl
911  << "bEquQvalue = " << bEquQvalue << std::endl
912  << "dbQ_db = " << dbQ_db << std::endl
913 
914  << "M_Inf = " << M_Inf<double>( v1Now ) << std::endl
915  << "M_Tau = " << M_Tau<double>( v1Now ) << std::endl
916  << "M_EquFvalue = " << M_EquFvalue << std::endl
917  << "dMF_dV1 = " << dMF_dV1 << std::endl
918  << "dMF_dM = " << dMF_dM << std::endl
919  << "M_EquQvalue = " << M_EquQvalue << std::endl
920  << "dMQ_dM = " << dMQ_dM << std::endl
921 
922  << "H_Inf = " << H_Inf<double>( v1Now ) << std::endl
923  << "H_Tau = " << H_Tau<double>( v1Now ) << std::endl
924  << "H_EquFvalue = " << H_EquFvalue << std::endl
925  << "dHF_dV1 = " << dHF_dV1 << std::endl
926  << "dHF_dH = " << dHF_dH << std::endl
927  << "H_EquQvalue = " << H_EquQvalue << std::endl
928  << "dHQ_dH = " << dHQ_dH << std::endl
929 
930  << "cEquFvalue = " << cEquFvalue << std::endl
931  << "dcF_dV1 = " << dcF_dV1 << std::endl
932  << "dcF_dc = " << dcF_dc << std::endl
933  << "cEquQvalue = " << cEquQvalue << std::endl
934  << "dcQ_dc = " << dcQ_dc << std::endl
935 
936  << "CaEquFvalue = " << CaEquFvalue << std::endl
937  << "dCaF_dV1 = " << dCaF_dV1 << std::endl
938  << "dCaF_dV2 = " << dCaF_dV2 << std::endl
939  << "dCaF_dM = " << dCaF_dM << std::endl
940  << "dCaF_dH = " << dCaF_dH << std::endl
941  << "dCaF_dCa = " << dCaF_dCa << std::endl
942  << "CaEquQvalue = " << CaEquQvalue << std::endl
943  << "dCaQ_dCa = " << dCaQ_dCa << std::endl
944 
945  << std::endl;
946  }
947 
948 
949  return bsuccess;
950 }
951 //-----------------------------------------------------------------------------
952 // Function : Instance::updatePrimaryState
953 // Purpose :
954 // Special Notes :
955 // Scope : public
956 // Creator : Richard Schiek, Electrical and Microsytem Modeling
957 // Creation Date : 01/02/08
958 //-----------------------------------------------------------------------------
960 {
961  bool bsuccess = true;
962 
964 
965  Linear::Vector * solVectorPtr = extData.nextSolVectorPtr;
966  Linear::Vector * staVectorPtr = extData.nextStaVectorPtr;
967 
968  (*staVectorPtr)[li_KCurrentState] = potassiumCurrent;
969  (*staVectorPtr)[li_NaCurrentState] = sodiumCurrent;
970 
971  return bsuccess;
972 }
973 
974 //-----------------------------------------------------------------------------
975 // Function : Instance::updateSecondaryState
976 // Purpose :
977 // Special Notes :
978 // Scope : public
979 // Creator : Richard Schiek, Electrical and Microsytem Modeling
980 // Creation Date : 01/02/08
981 //-----------------------------------------------------------------------------
983 {
984  bool bsuccess = true;
985 
986  return bsuccess;
987 }
988 
989 //-----------------------------------------------------------------------------
990 // Function : Instance::loadDAEQVector
991 //
992 // Purpose : Loads the Q-vector contributions for a single
993 // Neuron2 instance.
994 //
995 // Special Notes : The "Q" vector is part of a standard DAE formalism in
996 // which the system of equations is represented as:
997 //
998 // f(x) = dQ(x)/dt + F(x) - B(t) = 0
999 //
1000 // Scope : public
1001 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1002 // Creation Date : 01/02/08
1003 //-----------------------------------------------------------------------------
1005 {
1006  bool bsuccess = true;
1007 
1008  Linear::Vector * daeQVecPtr = extData.daeQVectorPtr;
1009 
1010  (*daeQVecPtr)[li_Pos] += kcl1Qvalue;
1011  (*daeQVecPtr)[li_Neg] += kcl2Qvalue;
1012  (*daeQVecPtr)[li_nPro] += nEquQvalue;
1013  (*daeQVecPtr)[li_mPro] += mEquQvalue;
1014  (*daeQVecPtr)[li_hPro] += hEquQvalue;
1015  (*daeQVecPtr)[li_aPro] += aEquQvalue;
1016  (*daeQVecPtr)[li_bPro] += bEquQvalue;
1017  (*daeQVecPtr)[li_M_Pro] += M_EquQvalue;
1018  (*daeQVecPtr)[li_H_Pro] += H_EquQvalue;
1019  (*daeQVecPtr)[li_cPro] += cEquQvalue;
1020  (*daeQVecPtr)[li_CaPro] += CaEquQvalue;
1021 
1022  return bsuccess;
1023 }
1024 
1025 
1026 
1027 //-----------------------------------------------------------------------------
1028 // Function : Instance::loadDAEFVector
1029 //
1030 // Purpose : Loads the F-vector contributions for a single
1031 // Neuron2 instance.
1032 //
1033 // Special Notes :
1034 //
1035 // Scope : public
1036 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1037 // Creation Date : 01/02/08
1038 //-----------------------------------------------------------------------------
1040 {
1041  bool bsuccess=true;
1042 
1043  Linear::Vector * daeFVecPtr = extData.daeFVectorPtr;
1044 
1045  (*daeFVecPtr)[li_Pos] += kcl1Fvalue;
1046  (*daeFVecPtr)[li_Neg] += kcl2Fvalue;
1047  (*daeFVecPtr)[li_nPro] += nEquFvalue;
1048  (*daeFVecPtr)[li_mPro] += mEquFvalue;
1049  (*daeFVecPtr)[li_hPro] += hEquFvalue;
1050  (*daeFVecPtr)[li_aPro] += aEquFvalue;
1051  (*daeFVecPtr)[li_bPro] += bEquFvalue;
1052  (*daeFVecPtr)[li_M_Pro] += M_EquFvalue;
1053  (*daeFVecPtr)[li_H_Pro] += H_EquFvalue;
1054  (*daeFVecPtr)[li_cPro] += cEquFvalue;
1055  (*daeFVecPtr)[li_CaPro] += CaEquFvalue;
1056 
1057  return bsuccess;
1058 }
1059 
1060 //-----------------------------------------------------------------------------
1061 // Function : Instance::loadDAEdQdx
1062 //
1063 // Purpose : Loads the Q-vector contributions for a single
1064 // Neuron 2 instance.
1065 // Scope : public
1066 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1067 // Creation Date : 01/02/08
1068 //-----------------------------------------------------------------------------
1070 {
1071  bool bsuccess = true;
1072 
1073  Linear::Matrix * dQdxMatPtr = extData.dQdxMatrixPtr;
1074 
1075  (*dQdxMatPtr)[li_Pos][APosEquPosNodeOffset] += dkcl1Q_dV1;
1076  (*dQdxMatPtr)[li_Pos][APosEquNegNodeOffset] += dkcl1Q_dV2;
1077 
1078  (*dQdxMatPtr)[li_Neg][ANegEquPosNodeOffset] += dkcl2Q_dV1;
1079  (*dQdxMatPtr)[li_Neg][ANegEquNegNodeOffset] += dkcl2Q_dV2;
1080 
1081  (*dQdxMatPtr)[li_nPro][ANEquNNodeOffset] += dnQ_dn;
1082  (*dQdxMatPtr)[li_mPro][AMEquMNodeOffset] += dmQ_dm;
1083  (*dQdxMatPtr)[li_hPro][AHEquHNodeOffset] += dhQ_dh;
1084  (*dQdxMatPtr)[li_aPro][AAEquANodeOffset] += daQ_da;
1085  (*dQdxMatPtr)[li_bPro][ABEquBNodeOffset] += dbQ_db;
1086  (*dQdxMatPtr)[li_M_Pro][AM_EquM_NodeOffset] += dMQ_dM;
1087  (*dQdxMatPtr)[li_H_Pro][AH_EquH_NodeOffset] += dHQ_dH;
1088  (*dQdxMatPtr)[li_cPro][ACEquCNodeOffset] += dcQ_dc;
1089  (*dQdxMatPtr)[li_CaPro][ACaEquCaNodeOffset] += dCaQ_dCa;
1090 
1091  return bsuccess;
1092 }
1093 
1094 
1095 
1096 //-----------------------------------------------------------------------------
1097 // Function : Instance::loadDAEdFdx ()
1098 //
1099 // Purpose : Loads the F-vector contributions for a single
1100 // Neuron 2 instance.
1101 //
1102 // Special Notes : This is an algebraic constaint.
1103 //
1104 // Scope : public
1105 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1106 // Creation Date : 01/02/08
1107 //-----------------------------------------------------------------------------
1109 {
1110  bool bsuccess = true;
1111 
1112  Linear::Matrix * dFdxMatPtr = extData.dFdxMatrixPtr;
1113 
1114  (*dFdxMatPtr)[li_Pos][APosEquPosNodeOffset] += dkcl1F_dV1;
1115  (*dFdxMatPtr)[li_Pos][APosEquNegNodeOffset] += dkcl1F_dV2;
1116  (*dFdxMatPtr)[li_Pos][APosEquNNodeOffset] += dkcl1F_dn;
1117  (*dFdxMatPtr)[li_Pos][APosEquMNodeOffset] += dkcl1F_dm;
1118  (*dFdxMatPtr)[li_Pos][APosEquHNodeOffset] += dkcl1F_dh;
1119  (*dFdxMatPtr)[li_Pos][APosEquANodeOffset] += dkcl1F_da;
1120  (*dFdxMatPtr)[li_Pos][APosEquBNodeOffset] += dkcl1F_db;
1121  (*dFdxMatPtr)[li_Pos][APosEquM_NodeOffset] += dkcl1F_dM;
1122  (*dFdxMatPtr)[li_Pos][APosEquH_NodeOffset] += dkcl1F_dH;
1123  (*dFdxMatPtr)[li_Pos][APosEquCNodeOffset] += dkcl1F_dc;
1124 
1125  (*dFdxMatPtr)[li_Neg][ANegEquPosNodeOffset] += dkcl2F_dV1;
1126  (*dFdxMatPtr)[li_Neg][ANegEquNegNodeOffset] += dkcl2F_dV2;
1127  (*dFdxMatPtr)[li_Neg][ANegEquNNodeOffset] += dkcl2F_dn;
1128  (*dFdxMatPtr)[li_Neg][ANegEquMNodeOffset] += dkcl2F_dm;
1129  (*dFdxMatPtr)[li_Neg][ANegEquHNodeOffset] += dkcl2F_dh;
1130  (*dFdxMatPtr)[li_Neg][ANegEquANodeOffset] += dkcl2F_da;
1131  (*dFdxMatPtr)[li_Neg][ANegEquBNodeOffset] += dkcl2F_db;
1132  (*dFdxMatPtr)[li_Neg][ANegEquM_NodeOffset] += dkcl2F_dM;
1133  (*dFdxMatPtr)[li_Neg][ANegEquH_NodeOffset] += dkcl2F_dH;
1134  (*dFdxMatPtr)[li_Neg][ANegEquCNodeOffset] += dkcl2F_dc;
1135 
1136  (*dFdxMatPtr)[li_nPro][ANEquPosNodeOffset] += dnF_dV1;
1137  (*dFdxMatPtr)[li_nPro][ANEquNNodeOffset] += dnF_dn;
1138 
1139  (*dFdxMatPtr)[li_mPro][AMEquPosNodeOffset] += dmF_dV1;
1140  (*dFdxMatPtr)[li_mPro][AMEquMNodeOffset] += dmF_dm;
1141 
1142  (*dFdxMatPtr)[li_hPro][AHEquPosNodeOffset] += dhF_dV1;
1143  (*dFdxMatPtr)[li_hPro][AHEquHNodeOffset] += dhF_dh;
1144 
1145  (*dFdxMatPtr)[li_aPro][AAEquPosNodeOffset] += daF_dV1;
1146  (*dFdxMatPtr)[li_aPro][AAEquANodeOffset] += daF_da;
1147 
1148  (*dFdxMatPtr)[li_bPro][ABEquPosNodeOffset] += dbF_dV1;
1149  (*dFdxMatPtr)[li_bPro][ABEquBNodeOffset] += dbF_db;
1150 
1151  (*dFdxMatPtr)[li_M_Pro][AM_EquPosNodeOffset] += dMF_dV1;
1152  (*dFdxMatPtr)[li_M_Pro][AM_EquM_NodeOffset] += dMF_dM;
1153 
1154  (*dFdxMatPtr)[li_H_Pro][AH_EquPosNodeOffset] += dHF_dV1;
1155  (*dFdxMatPtr)[li_H_Pro][AH_EquH_NodeOffset] += dHF_dH;
1156 
1157  (*dFdxMatPtr)[li_cPro][ACEquPosNodeOffset] += dcF_dV1;
1158  (*dFdxMatPtr)[li_cPro][ACEquCNodeOffset] += dcF_dc;
1159  (*dFdxMatPtr)[li_cPro][ACEquCaNodeOffset] += dcF_dCa;
1160 
1161  (*dFdxMatPtr)[li_CaPro][ACaEquPosNodeOffset] += dCaF_dV1;
1162  (*dFdxMatPtr)[li_CaPro][ACaEquNegNodeOffset] += dCaF_dV2;
1163  (*dFdxMatPtr)[li_CaPro][ACaEquM_NodeOffset] += dCaF_dM;
1164  (*dFdxMatPtr)[li_CaPro][ACaEquH_NodeOffset] += dCaF_dH;
1165  (*dFdxMatPtr)[li_CaPro][ACaEquCaNodeOffset] += dCaF_dCa;
1166 
1167  return bsuccess;
1168 }
1169 
1170 //-----------------------------------------------------------------------------
1171 // Function : Instance::setIC
1172 // Purpose :
1173 // Special Notes :
1174 // Scope : public
1175 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1176 // Creation Date : 01/02/08
1177 //-----------------------------------------------------------------------------
1179 {
1180  bool bsuccess = true;
1181 
1182  return bsuccess;
1183 }
1184 
1185 //-----------------------------------------------------------------------------
1186 // Function : Instance::varTypes
1187 // Purpose :
1188 // Special Notes :
1189 // Scope : public
1190 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1191 // Creation Date : 01/02/08
1192 //-----------------------------------------------------------------------------
1193 void Instance::varTypes( std::vector<char> & varTypeVec )
1194 {
1195  //varTypeVec.resize(1);
1196  //varTypeVec[0] = 'I';
1197 }
1198 
1199 
1200 //-----------------------------------------------------------------------------
1201 // Function : Model::processParams
1202 // Purpose :
1203 // Special Notes :
1204 // Scope : public
1205 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1206 // Creation Date : 01/02/08
1207 //-----------------------------------------------------------------------------
1209 {
1210  return true;
1211 }
1212 
1213 //----------------------------------------------------------------------------
1214 // Function : Model::processInstanceParams
1215 // Purpose :
1216 // Special Notes :
1217 // Scope : public
1218 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1219 // Creation Date : 01/02/08
1220 //----------------------------------------------------------------------------
1222 {
1223 
1224  std::vector<Instance*>::iterator iter;
1225  std::vector<Instance*>::iterator first = instanceContainer.begin();
1226  std::vector<Instance*>::iterator last = instanceContainer.end();
1227 
1228  for (iter=first; iter!=last; ++iter)
1229  {
1230  (*iter)->processParams();
1231  }
1232 
1233  return true;
1234 }
1235 
1236 //-----------------------------------------------------------------------------
1237 // Function : Model::Model
1238 // Purpose : block constructor
1239 // Special Notes :
1240 // Scope : public
1241 // Creator : Richard Schiek, Electrical and Microsytem Modeling
1242 // Creation Date : 01/02/08
1243 //-----------------------------------------------------------------------------
1245  const Configuration & configuration,
1246  const ModelBlock & MB,
1247  const FactoryBlock & factory_block)
1248  : DeviceModel(MB, configuration.getModelParameters(), factory_block)
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 : 01/02/08
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 : 01/02/08
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", 2)
1354  .registerModelType("neuron", 2);
1355 }
1356 
1357 } // namespace Neuron2
1358 } // namespace Device
1359 } // namespace Xyce
const InstanceName & getName() const
void registerLIDs(const std::vector< int > &intLIDVecRef, const std::vector< int > &extLIDVecRef)
static ScalarT C_EquF(const ScalarT &Vn1, const ScalarT &C, const ScalarT &CaConc, const ScalarT &Vrest)
static ScalarT H_EquF(const ScalarT &Vn1, const ScalarT &H, const ScalarT &Vrest)
const SolverState & solverState_
Descriptor & addPar(const char *parName, T default_value, T U::*varPtr)
Adds the parameter description to the parameter map.
Definition: N_DEV_Pars.h:1429
static Device * factory(const Configuration &configuration, const FactoryBlock &factory_block)
void registerStateLIDs(const std::vector< int > &staLIDVecRef)
Linear::Vector * nextSolVectorPtr
Linear::Vector * daeQVectorPtr
Pure virtual class to augment a linear system.
void addInternalNode(Util::SymbolTable &symbol_table, int index, const InstanceName &instance_name, const std::string &lead_name)
static ScalarT H_EquQ(const ScalarT &H)
static ScalarT kcl2EquF(const ScalarT &Vn1, const ScalarT &Vn2, const ScalarT &n, const ScalarT &m, const ScalarT &h, const ScalarT &a, const ScalarT &b, const ScalarT &MC, const ScalarT &HC, const ScalarT &CC, const ScalarT &memG, const ScalarT &restV, const ScalarT &Kg, const ScalarT &Ke, const ScalarT &NaG, const ScalarT &NaE, const ScalarT &Ag, const ScalarT &Ae, const ScalarT &CaTg, const ScalarT &CaE, const ScalarT &KCaG)
static void loadModelParameters(ParametricData< Model > &model_parameters)
Definition: N_DEV_Neuron2.C:75
#define AssertLIDs(cmp)
static ScalarT mEquF(const ScalarT &Vn1, const ScalarT &m, const ScalarT &Vrest)
bool processInstanceParams()
processInstanceParams
virtual void registerJacLIDs(const JacobianStamp &jacLIDVec)
static ScalarT nEquQ(const ScalarT &n)
static ScalarT hEquQ(const ScalarT &h)
void loadNodeSymbols(Util::SymbolTable &symbol_table) const
Populates and returns the store name map.
static ScalarT aEquF(const ScalarT &Vn1, const ScalarT &a, const ScalarT &Vrest)
std::vector< Param > params
Parameters from the line.
static void loadInstanceParameters(ParametricData< Instance > &instance_parameters)
Definition: N_DEV_Neuron2.C:71
const std::string & getName() const
static ScalarT aEquQ(const ScalarT &a)
bool processParams()
processParams
static ScalarT Ca_EquQ(const ScalarT &Ca)
The FactoryBlock contains parameters needed by the device, instance and model creation functions...
const DeviceOptions & deviceOptions_
static ScalarT kcl1EquQ(const ScalarT &Vn1, const ScalarT &Vn2, const ScalarT &memC)
static ScalarT kcl1EquF(const ScalarT &Vn1, const ScalarT &Vn2, const ScalarT &n, const ScalarT &m, const ScalarT &h, const ScalarT &a, const ScalarT &b, const ScalarT &MC, const ScalarT &HC, const ScalarT &CC, const ScalarT &memG, const ScalarT &restV, const ScalarT &Kg, const ScalarT &Ke, const ScalarT &NaG, const ScalarT &NaE, const ScalarT &Ag, const ScalarT &Ae, const ScalarT &CaTg, const ScalarT &CaE, const ScalarT &KCaG)
static ScalarT M_EquF(const ScalarT &Vn1, const ScalarT &M, const ScalarT &Vrest)
bool updateTemperature(const double &temp_tmp)
Linear::Vector * nextStaVectorPtr
static Config< T > & addConfiguration()
Adds the device to the Xyce device configuration.
Linear::Matrix * dFdxMatrixPtr
static ScalarT bEquF(const ScalarT &Vn1, const ScalarT &b, const ScalarT &Vrest)
virtual void forEachInstance(DeviceInstanceOp &op) const
Apply a device instance "op" to all instances associated with this model.
The Device class is an interface for device implementations.
Definition: N_DEV_Device.h:101
static ScalarT nEquF(const ScalarT &Vn1, const ScalarT &n, const ScalarT &Vrest)
Class Configuration contains device configuration data.
Instance(const Configuration &configuration, const InstanceBlock &IB, Model &Miter, const FactoryBlock &factory_block)
const SolverState & getSolverState() const
static ScalarT Ca_EquF(const ScalarT &Vn1, const ScalarT &Vn2, const ScalarT &MC, const ScalarT &HC, const ScalarT &Ca, const ScalarT &CaTg, const ScalarT &CaE, const ScalarT &CaGamma, const ScalarT &CaTau)
static ScalarT C_EquQ(const ScalarT &C)
void registerJacLIDs(const std::vector< std::vector< int > > &jacLIDVec)
Linear::Vector * daeFVectorPtr
static ScalarT mEquQ(const ScalarT &m)
std::vector< Instance * > instanceContainer
static std::vector< std::vector< int > > jacStamp
Definition: N_DEV_Neuron2.h:95
virtual std::ostream & printOutInstances(std::ostream &os) const
ModelBlock represents a .MODEL line from the netlist.
Manages parameter binding for class C.
Definition: N_DEV_Pars.h:214
InstanceBlock represent a device instance line from the netlist.
Linear::Matrix * dQdxMatrixPtr
const std::vector< std::vector< int > > & jacobianStamp() const
void varTypes(std::vector< char > &varTypeVec)
static ScalarT kcl2EquQ(const ScalarT &Vn1, const ScalarT &Vn2, const ScalarT &memC)
void setModParams(const std::vector< Param > &params)
static ScalarT hEquF(const ScalarT &Vn1, const ScalarT &h, const ScalarT &Vrest)