Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_DEV_Neuron3.h
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$
27 //
28 // Purpose : Neuron classes.
29 //
30 // Special Notes :
31 //
32 // Creator : Richard Schiek, SNL, Electrical and Microsystem Modeling
33 //
34 // Creation Date : 06/10/09
35 //
36 // Revision Information:
37 // ---------------------
38 //
39 // Revision Number: $Revision$
40 //
41 // Revision Date : $Date$
42 //
43 // Current Owner : $Author$
44 //-----------------------------------------------------------------------------
45 
46 #ifndef Xyce_N_DEV_Neuron3_h
47 #define Xyce_N_DEV_Neuron3_h
48 
49 #include <N_DEV_Configuration.h>
50 #include <Sacado.hpp>
51 
52 // ---------- Xyce Includes ----------
53 #include <N_DEV_DeviceBlock.h>
54 #include <N_DEV_DeviceInstance.h>
55 #include <N_DEV_DeviceModel.h>
56 
57 #include <N_DEV_Neuron.h>
58 
59 #ifdef HAVE_MATH_H
60 #include <math.h>
61 #endif
62 
63 namespace Xyce {
64 namespace Device {
65 namespace Neuron3 {
66 
67 class Model;
68 class Instance;
69 
70 struct Traits : public DeviceTraits<Model, Instance, Neuron::Traits>
71 {
72  static const char *name() {return "Neuron";}
73  static const char *deviceTypeName() {return "YNEURON level 3";}
74  static const int numNodes() {return 2;}
75  static const bool modelRequired() {return true;}
76  static const bool isLinearDevice() {return true;}
77 
78  static Device *factory(const Configuration &configuration, const FactoryBlock &factory_block);
79  static void loadModelParameters(ParametricData<Model> &model_parameters);
80  static void loadInstanceParameters(ParametricData<Instance> &instance_parameters);
81 };
82 
83 //-----------------------------------------------------------------------------
84 // Class : Instance
85 // Purpose : This is class refers to a single instance of the
86 // Neuron device. It has two nodes associated with it, a
87 // positive and a negative node. See the NeuronInstance
88 // class for a more detailed explanation.
89 // Special Notes :
90 // Creator : Richard Schiek, SNL, Electrical and Microsystem Modeling
91 // Creation Date : 06/10/09
92 //-----------------------------------------------------------------------------
93 class Instance : public DeviceInstance
94 {
95  friend class ParametricData<Instance>;
96  friend class Model;
97  friend class Traits;
98 
99 public:
100 
101  Instance(
102  const Configuration & configuration,
103  const InstanceBlock & IB,
104  Model & Miter,
105  const FactoryBlock & factory_block);
106 
107 
108  ~Instance();
109 
110 private:
111  Instance(const Instance &);
112  Instance &operator=(const Instance &);
113 
114 public:
115  void registerLIDs( const std::vector<int> & intLIDVecRef,
116  const std::vector<int> & extLIDVecRef );
117  void registerStateLIDs( const std::vector<int> & staLIDVecRef );
118 
119  std::map<int,std::string> & getIntNameMap ();
120  bool loadDeviceMask();
121  const std::vector< std::vector<int> > & jacobianStamp() const;
122  void registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec );
123 
124  bool processParams ();
125  bool updateTemperature(const double & temp_tmp);
126 
127  bool updateIntermediateVars ();
128  bool updatePrimaryState ();
129  bool updateSecondaryState ();
130  bool setIC ();
131 
132  void varTypes( std::vector<char> & varTypeVec );
133 
134  // load functions, residual:
135  bool loadDAEQVector ();
136  bool loadDAEFVector ();
137 
138  void auxDAECalculations ();
139 
140  // load functions, Jacobian:
141  bool loadDAEdQdx ();
142  bool loadDAEdFdx ();
143 
144 protected:
145 private:
146  // These functions represents the equations that need to be solved
147  // for this device. Since Xyce loads an F and Q contribution, the
148  // equations are broken up into their F and Q components. Thus there
149  // is a kcl1EquF() and kcl1EquQ(). Automatic differentiation will
150  // be used to generate all the derivatives of these equations for the
151  // dF/dX and dQ/dX loads
152 
153  // first we list some utility functions for calculating coefficients
154  //
155  // These functions expect V to be in milli-volts and then return values that
156  // are in 1/ms. Thus the extra factor's of 1000 here and there
157 
158  // potassium current, functions for activator equation
159  template <typename ScalarT>
160  static ScalarT alphaN( const ScalarT & Vn1, const ScalarT & Vrest)
161  {
162  ScalarT vDiff = 1000.0 * (Vn1 - Vrest); // convert voltage to milli-volts
163  ScalarT r;
164  // if ((vDiff > 9.99) && (vDiff < 10.01) )
165  // {
166  // r = 1.0/(10.0 * ( std::exp( (10.0 - vDiff)/10.0 )));
167  // }
168  // else
169  // {
170  r = (10.0 - vDiff) /
171  (100.0 * ( std::exp( (10.0 - vDiff)/10.0 ) - 1.0 ));
172  // }
173  r *= 1000.0; // change from 1/ms to 1/s
174  return r;
175  }
176 
177  template <typename ScalarT>
178  static ScalarT betaN( const ScalarT & Vn1, const ScalarT & Vrest)
179  {
180  ScalarT vDiff = 1000.0 * (Vn1 - Vrest);
181  ScalarT r = 0.125 * std::exp( -vDiff/80.0 );
182  r *= 1000.0; // change from 1/ms to 1/s
183  return r;
184  }
185 
186  // sodium current, functions for activator equation
187  template <typename ScalarT>
188  static ScalarT alphaM( const ScalarT & Vn1, const ScalarT & Vrest)
189  {
190  ScalarT vDiff = 1000.0 * (Vn1 - Vrest);
191  ScalarT r;
192  // if ((vDiff > 24.99) && (vDiff < 25.01) )
193  // {
194  // r = (1.0) /
195  // (( std::exp( (25.0 - vDiff)/10.0 )));
196  // }
197  // else
198  // {
199  r = (25.0 - vDiff) /
200  (10.0 * ( std::exp( (25.0 - vDiff)/10.0 ) - 1.0 ));
201  // }
202  r *= 1000.0; // change from 1/ms to 1/s
203  return r;
204  }
205 
206  template <typename ScalarT>
207  static ScalarT betaM( const ScalarT & Vn1, const ScalarT & Vrest)
208  {
209  ScalarT vDiff = 1000.0 * (Vn1 - Vrest);
210  ScalarT r = 4.0 * std::exp( -vDiff/18.0 );
211  r *= 1000.0; // change from 1/ms to 1/s
212 
213  return r;
214  }
215 
216  template <typename ScalarT>
217  static ScalarT alphaH( const ScalarT & Vn1, const ScalarT & Vrest)
218  {
219  ScalarT vDiff = 1000.0 * (Vn1 - Vrest);
220  ScalarT r = 0.07 * std::exp( -vDiff/20.0 );
221  r *= 1000.0; // change from 1/ms to 1/s
222  return r;
223  }
224 
225  template <typename ScalarT>
226  static ScalarT betaH( const ScalarT & Vn1, const ScalarT & Vrest)
227  {
228  ScalarT vDiff = 1000.0 * (Vn1 - Vrest);
229  ScalarT r = 1.0 / ( std::exp( (30.0 - vDiff)/10.0 ) + 1.0 );
230  r *= 1000.0; // change from 1/ms to 1/s
231  return r;
232  }
233 
234  // now the device equations
235  // KCL equation 1
236  // VSeg -- segment voltage
237  // VSegP -- previous segment voltage
238  // VSegN -- next segment voltage
239  // n, m, h -- vars which determine current through themembrane
240  // gPrev, gNext -- conductivity to previous and next segment (can be zero)
241  //
242  // Full equation is:
243  // i_trans_membrane - i_externally_injected - gNext (VSegN - VSeg) - gPrev (VSegP - VSeg ) + Cmem dVseg/dt = 0
244  //
245  // i_tarns_membrane is determined by the ion channel equations -- hodgkin-huxley in this case.
246  // i_externally_injeted is not supported here but could be an externally injected current.
247  //
248  template <typename ScalarT>
249  static ScalarT kcl1EquF( const ScalarT& VSeg, const ScalarT& VSegP, const ScalarT & VSegN, const ScalarT& n, const ScalarT& m, const ScalarT& h,
250  const ScalarT& gPrev, const ScalarT& gNext, const ScalarT& memG, const ScalarT& restV, const ScalarT& Kg, const ScalarT& Ke, const ScalarT& NaG, const ScalarT& NaE )
251  {
252  ScalarT powN = n * n * n * n;
253  ScalarT powM = m * m * m;
254  //ScalarT r = memG * (VSeg - restV) + Kg * powN * (VSeg - Ke ) + NaG * powM * h * (VSeg - NaE ) - gNext * (VSegN - VSeg) - gPrev * (VSegP - VSeg);
255  ScalarT r = memG * (VSeg - restV) + Kg * powN * (VSeg - Ke ) + NaG * powM * h * (VSeg - NaE )
256  - gNext * (VSegN - VSeg) - gPrev * (VSegP - VSeg);
257  return r;
258  }
259 
260  template <typename ScalarT>
261  static ScalarT kcl1EquQ( const ScalarT& VSeg, const ScalarT& memC )
262  {
263  ScalarT r = memC * VSeg;
264  return r;
265  }
266 
267 #if 0
268  // KCL equation 2 -- -1 * equation 1 because of device symmetry
269  template <typename ScalarT>
270  static ScalarT kcl2EquF( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& n, const ScalarT& m, const ScalarT& h,
271  const ScalarT& memG, const ScalarT& restV, const ScalarT& Kg, const ScalarT& Ke, const ScalarT& NaG, const ScalarT& NaE )
272  {
273  ScalarT powN = n * n * n * n;
274  ScalarT powM = m * m * m;
275  ScalarT r = -1.0*(memG * (Vn1 - Vn2 - restV) + Kg * powN * (Vn1 - Vn2 - Ke ) + NaG * powM * h * (Vn1 - Vn2 - NaE ));
276  //ScalarT r = -1.0*(memG * (Vn1 - restV) - Kg * powN * (Vn1 - Vn2 - Ke ) - NaG * powM * h * (Vn1 - Vn2 - NaE ));
277  return r;
278  }
279 
280  template <typename ScalarT>
281  static ScalarT kcl2EquQ( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& memC )
282  {
283  ScalarT r = -1.0 * memC * (Vn1 - Vn2);
284  return r;
285  }
286 #endif
287 
288  // n conservation equation
289  template <typename ScalarT>
290  static ScalarT nEquF( const ScalarT& Vn1, const ScalarT& n, const ScalarT& Vrest )
291  {
292  ScalarT r = alphaN<ScalarT>( Vn1, Vrest ) * (1.0 - n ) - betaN<ScalarT>( Vn1, Vrest ) * n;
293  return r;
294  }
295 
296  template <typename ScalarT>
297  static ScalarT nEquQ( const ScalarT& n )
298  {
299  ScalarT r = -n;
300  return r;
301  }
302 
303  // m conservation equation
304  template <typename ScalarT>
305  static ScalarT mEquF( const ScalarT& Vn1, const ScalarT& m, const ScalarT& Vrest )
306  {
307  ScalarT r = alphaM<ScalarT>( Vn1, Vrest ) * (1.0 - m ) - betaM<ScalarT>( Vn1, Vrest ) * m;
308  return r;
309  }
310 
311  template <typename ScalarT>
312  static ScalarT mEquQ( const ScalarT& m )
313  {
314  ScalarT r = -m;
315  return r;
316  }
317 
318  // h conservation equation
319  template <typename ScalarT>
320  static ScalarT hEquF( const ScalarT& Vn1, const ScalarT& h, const ScalarT& Vrest )
321  {
322  ScalarT r = alphaH<ScalarT>( Vn1, Vrest ) * (1.0 - h ) - betaH<ScalarT>( Vn1, Vrest ) * h;
323  return r;
324  }
325 
326  template <typename ScalarT>
327  static ScalarT hEquQ( const ScalarT& h )
328  {
329  ScalarT r = -h;
330  return r;
331  }
332 
333 public:
334  // iterator reference to the Neuron model which owns this instance.
335  // Getters and setters
337  {
338  return model_;
339  }
340 
341 private:
342 
343  Model & model_; //< Owning model
344 
345  std::vector< std::vector<int> > jacStamp;
346 
347  // model level parameters that can be overridden at the instance level
348  double rInt; // intracellular resistivity
349  double radius; // Segment radius
350  double length; // cable length (segment length = length/nSeg)
351  double segArea; // segment area (derrived from radius, length and nSeg)
352  int nSeg; // number of segments
353  bool rIntGiven;
356  bool nSegGiven;
357 
358  // the cable equation is dependent on parameters from the previous and next segment of other
359  // cables attached to the input/ouput nodes (i.e. when we branch a cable). This lets
360  // one set the parameters for the previous/next cable segment.
361  double rIntPrevious;
364  double rIntNext;
365  double radiusNext;
366  double lengthNext;
373 
374  // conductance between segments -- calculated from radius, rInt and length and number of segments
375  // gBackward connects to previous node. gForward connects to the next node
376  std::vector<double> gBackward;
377  std::vector<double> gForward;
378 
379  // derrived quantities computed in updateIntermediateVars
380  // and used in the load functions (no q terms on the external nodes)
381  double kcl1Fvalue;
382  double kcl2Fvalue;
383  // internal segments
384  std::vector<double> segFvalue;
385  std::vector<double> segQvalue;
386  std::vector<double> segNEquFvalue, segNEquQvalue;
387  std::vector<double> segMEquFvalue, segMEquQvalue;
388  std::vector<double> segHEquFvalue, segHEquQvalue;
389  // jacobian terms
392  // internal equations
393  std::vector<double> segF_dVp, segF_dV, segF_dVn, segF_dn, segF_dm, segF_dh;
394  std::vector<double> segQ_dV;
395  std::vector<double> dnF_dV, dnF_dn, dnQ_dn;
396  std::vector<double> dmF_dV, dmF_dm, dmQ_dm;
397  std::vector<double> dhF_dV, dhF_dh, dhQ_dh;
398 
399  // state variables
400  std::vector<double> potassiumCurrent;
401  std::vector<double> sodiumCurrent;
402 
403  // local state indices (offsets)
404  std::vector<int> li_KCurrentState;
405  std::vector<int> li_NaCurrentState;
406 
407  // local solution indices (offsets)
408  int li_Pos; // local index to positive node on this device
409  int li_Neg; // local index to negative node on this device
410  // local solution indices for internal vars (variable number of these)
411  std::vector<int> li_Vol; // local index to segment voltage
412  std::vector<int> li_nPro; // local index to n promoter value (Na current)
413  std::vector<int> li_mPro; // local index to m promoter value (K current)
414  std::vector<int> li_hPro; // local index to h promoter value (K current)
415 
416  // Matrix equation index variables:
417 
418  // Offset variables corresponding to the above declared indices.
421  std::vector<int> SegVEqnVpreOffset;
422  std::vector<int> SegVEqnVsegOffset;
423  std::vector<int> SegVEqnVnexOffset;
424  std::vector<int> SegVEqnNOffset;
425  std::vector<int> SegVEqnMOffset;
426  std::vector<int> SegVEqnHOffset;
427  std::vector<int> NEquVNodeOffset;
428  std::vector<int> NEquNNodeOffset;
429  std::vector<int> MEquVNodeOffset;
430  std::vector<int> MEquMNodeOffset;
431  std::vector<int> HEquVNodeOffset;
432  std::vector<int> HEquHNodeOffset;
433 };
434 
435 //-----------------------------------------------------------------------------
436 // Class : Model
437 // Purpose :
438 // Special Notes :
439 // Creator : Richard Schiek, SNL, Electrical and Microsystem Modeling
440 // Creation Date : 06/10/09
441 //-----------------------------------------------------------------------------
442 class Model : public DeviceModel
443 {
444  typedef std::vector<Instance *> InstanceVector;
445 
446  friend class ParametricData<Model>;
447  friend class Instance;
448  friend class Traits;
449 
450 public:
451  Model(
452  const Configuration & configuration,
453  const ModelBlock & MB,
454  const FactoryBlock & factory_block);
455  ~Model();
456 
457 private:
458  Model();
459  Model(const Model &);
460  Model &operator=(const Model &);
461 
462 public:
463  virtual void forEachInstance(DeviceInstanceOp &op) const /* override */;
464 
465  virtual std::ostream &printOutInstances(std::ostream &os) const;
466 
467  bool processParams ();
468  bool processInstanceParams ();
469 
470 private:
471 
472  // parameter variables
473  double cMem; // membrane capacitance
474  double gMem; // membrane conductance
475  double vRest; // resting potential
476  double eNa; // sodium rest potential
477  double gNa; // sodium base conductance
478  double eK; // potassium rest potential
479  double gK; // potassium base conductance
480  double rInt; // intracellular resistivity
481  double radius; // Segment radius
482  double length; // cable length (segment length = length/nSeg)
483  int nSeg; // number of segments
484 
485  // the cable equation is dependent on parameters from the previous and next segment of other
486  // cables attached to the input/ouput nodes (i.e. when we branch a cable). This lets
487  // one set the parameters for the previous/next cable segment.
488  double rIntPrevious;
491  double rIntNext;
492  double radiusNext;
493  double lengthNext;
500 
501  // flags that parameters were given
502  bool cMemGiven;
503  bool gMemGiven;
505  bool eNaGiven;
506  bool gNaGiven;
507  bool eKGiven;
508  bool gKGiven;
509  bool rIntGiven;
512  bool nSegGiven;
513 
514 
515 public:
516  void addInstance(Instance *instance)
517  {
518  instanceContainer.push_back(instance);
519  }
520 
522  {
523  return instanceContainer;
524  }
525 
527  {
528  return instanceContainer;
529  }
530 
531 private:
532  std::vector<Instance*> instanceContainer;
533 };
534 
535 void registerDevice();
536 
537 } // namespace Neuron3
538 } // namespace Device
539 } // namespace Xyce
540 
543 
544 #endif