Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
N_DEV_Neuron9.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: N_DEV_Neuron9.h,v $
27 //
28 // Purpose : Neuron classes.
29 //
30 // Special Notes :
31 //
32 // Creator : Christy Warrender, SNL, Cognitive Modeling
33 //
34 // Creation Date : 06/22/12
35 //
36 // Revision Information:
37 // ---------------------
38 //
39 // Revision Number: $Revision: 1.21.2.1 $
40 //
41 // Revision Date : $Date: 2014/02/26 20:16:30 $
42 //
43 // Current Owner : $Author: tvrusso $
44 //-----------------------------------------------------------------------------
45 
46 #ifndef Xyce_N_DEV_Neuron9_h
47 #define Xyce_N_DEV_Neuron9_h
48 
49 #include <N_DEV_Configuration.h>
50 #include <Sacado.hpp>
51 
52 // ---------- Xyce Includes ----------
53 #include <N_DEV_DeviceMaster.h>
54 #include <N_DEV_DeviceBlock.h>
55 #include <N_DEV_DeviceInstance.h>
56 #include <N_DEV_DeviceModel.h>
57 
58 #include <N_DEV_Neuron.h>
59 
60 #ifdef HAVE_MATH_H
61 #include <math.h>
62 #endif
63 
64 namespace Xyce {
65 namespace Device {
66 namespace Neuron9 {
67 
68 class Model;
69 class Instance;
70 
71 struct Traits : public DeviceTraits<Model, Instance, Neuron::Traits>
72 {
73  static const char *name() {return "Neuron";}
74  static const char *deviceTypeName() {return "YNEURON level 9";}
75  static const int numNodes() {return 2;}
76  static const bool modelRequired() {return true;}
77  static const bool isLinearDevice() {return true;}
78 
79  static Device *factory(const Configuration &configuration, const FactoryBlock &factory_block);
80  static void loadModelParameters(ParametricData<Model> &model_parameters);
81  static void loadInstanceParameters(ParametricData<Instance> &instance_parameters);
82 };
83 
84 //-----------------------------------------------------------------------------
85 // Class : Instance
86 // Purpose : This is class refers to a single instance of the
87 // Neuron device. It has two nodes associated with it, a
88 // positive and a negative node. See the NeuronInstance
89 // class for a more detailed explanation.
90 // Special Notes :
91 // Creator : Christy Warrender, SNL, Cognitive Modeling
92 // Creation Date : 06/22/12
93 //-----------------------------------------------------------------------------
94 class Instance : public DeviceInstance
95 {
96  friend class ParametricData<Instance>;
97  friend class Model;
98  friend class Traits;friend class Master;
99 
100 public:
101  static std::vector< std::vector<int> > jacStamp;
102 
103  Instance(
104  const Configuration & configuration,
105  const InstanceBlock & IB,
106  Model & Miter,
107  const FactoryBlock & factory_block);
108 
109 
110  ~Instance();
111 
112 private:
113  Instance(const Instance &);
114  Instance &operator=(const Instance &);
115 
116 public:
117  void registerLIDs( const std::vector<int> & intLIDVecRef,
118  const std::vector<int> & extLIDVecRef );
119  void registerStateLIDs( const std::vector<int> & staLIDVecRef );
120 
121  std::map<int,std::string> & getIntNameMap ();
122  bool loadDeviceMask();
123  const std::vector< std::vector<int> > & jacobianStamp() const;
124  void registerJacLIDs( const std::vector< std::vector<int> > & jacLIDVec );
125 
126  bool processParams ();
127  bool updateTemperature(const double & temp_tmp);
128 
129  bool updateIntermediateVars ();
130  bool updatePrimaryState ();
131  bool updateSecondaryState ();
132  bool setIC ();
133 
134  void varTypes( std::vector<char> & varTypeVec );
135 
136  // load functions, residual:
137  bool loadDAEQVector ();
138  bool loadDAEFVector ();
139 
140  void auxDAECalculations ();
141 
142  // load functions, Jacobian:
143  bool loadDAEdQdx ();
144  bool loadDAEdFdx ();
145 
146 private:
147  // These functions represent the equations that need to be solved
148  // for this device. Since Xyce loads an F and Q contribution, the
149  // equations are broken up into their F and Q components. Thus there
150  // is a kcl1EquF() and kcl1EquQ(). Automatic differentiation will
151  // be used to generate all the derivatives of these equations for the
152  // dF/dX and dQ/dX loads
153 
154  // first we list some utility functions for calculating coefficients.
155  // alpha and beta equations are taken from Brette et al 07.
156  // They're generally functions of membrane voltage; here, the membrane voltage is
157  // the difference between Vn1 and Vn2.
158  // These functions expect V to be in milli-volts and then return values that
159  // are in 1/ms. Thus the extra factor's of 1000 here and there
160 
161  // potassium current, functions for activator equation
162  template <typename ScalarT>
163  static ScalarT alphaN( const ScalarT & Vn1, const ScalarT & Vn2, const ScalarT & Vrest)
164  {
165  ScalarT vDiff = 1000.0 * (Vn1 - Vn2); // convert voltage to milli-volts
166  ScalarT r;
167  r = 0.032*(15.0 - vDiff + VT)/( std::exp( (15.0 - vDiff + VT)/5.0) - 1.0);
168  r *= 1000.0; // change from 1/ms to 1/s
169  return r;
170  }
171 
172  template <typename ScalarT>
173  static ScalarT betaN( const ScalarT & Vn1, const ScalarT & Vn2, const ScalarT & Vrest)
174  {
175  ScalarT vDiff = 1000.0 * (Vn1 - Vn2); // convert voltage to milli-volts
176  ScalarT r;
177  r = 0.5*( std::exp( (10.0 - vDiff + VT)/40.0) );
178  r *= 1000.0; // change from 1/ms to 1/s
179  return r;
180  }
181 
182  // sodium current, functions for activator equation
183  template <typename ScalarT>
184  static ScalarT alphaM( const ScalarT & Vn1, const ScalarT & Vn2, const ScalarT & Vrest)
185  {
186  ScalarT vDiff = 1000.0 * (Vn1 - Vn2); // convert voltage to milli-volts
187  ScalarT r;
188  r = 0.32*(13.0 - vDiff + VT )/( std::exp((13.0 - vDiff + VT )/4.0) - 1.0);
189  r *= 1000.0; // change from 1/ms to 1/s
190  return r;
191  }
192 
193  template <typename ScalarT>
194  static ScalarT betaM( const ScalarT & Vn1, const ScalarT & Vn2, const ScalarT & Vrest)
195  {
196  ScalarT vDiff = 1000.0 * (Vn1 - Vn2); // convert voltage to milli-volts
197  ScalarT r;
198  r = 0.28*(vDiff - VT - 40.0)/( std::exp((vDiff - VT - 40.0)/5.0) - 1.0);
199  r *= 1000.0; // change from 1/ms to 1/s
200  return r;
201  }
202 
203  template <typename ScalarT>
204  static ScalarT alphaH( const ScalarT & Vn1, const ScalarT & Vn2, const ScalarT & Vrest)
205  {
206  ScalarT vDiff = 1000.0 * (Vn1 - Vn2); // convert voltage to milli-volts
207  ScalarT r;
208  r = 0.128 * std::exp( (17.0 - vDiff + VT)/18.0 );
209  r *= 1000.0; // change from 1/ms to 1/s
210  return r;
211  }
212 
213  template <typename ScalarT>
214  static ScalarT betaH( const ScalarT & Vn1, const ScalarT & Vn2, const ScalarT & Vrest)
215  {
216  ScalarT vDiff = 1000.0 * (Vn1 - Vn2); // convert voltage to milli-volts
217  ScalarT r;
218  r = 4.0/( 1.0 + std::exp( (40.0 - vDiff + VT)/5.0) );
219  r *= 1000.0; // change from 1/ms to 1/s
220  return r;
221  }
222 
223  // now the device equations
224  // KCL equation 1
225  template <typename ScalarT>
226  static ScalarT kcl1EquF( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& n, const ScalarT& m, const ScalarT& h,
227  const ScalarT& memG, const ScalarT& leakE, const ScalarT& Kg, const ScalarT& Ke, const ScalarT& NaG, const ScalarT& NaE )
228  {
229  ScalarT powN = n * n * n * n;
230  ScalarT powM = m * m * m;
231  ScalarT r = memG * (Vn1 - Vn2 - leakE) + Kg * powN * (Vn1 - Vn2 - Ke ) + NaG * powM * h * (Vn1 - Vn2 - NaE );
232  return r;
233  }
234 
235  template <typename ScalarT>
236  static ScalarT kcl1EquQ( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& memC )
237  {
238  ScalarT r = memC * (Vn1 - Vn2);
239  return r;
240  }
241 
242  // KCL equation 2 -- -1 * equation 1 because of device symmetry
243  template <typename ScalarT>
244  static ScalarT kcl2EquF( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& n, const ScalarT& m, const ScalarT& h,
245  const ScalarT& memG, const ScalarT& leakE, const ScalarT& Kg, const ScalarT& Ke, const ScalarT& NaG, const ScalarT& NaE )
246  {
247  ScalarT powN = n * n * n * n;
248  ScalarT powM = m * m * m;
249  ScalarT r = -1.0*(memG * (Vn1 - Vn2 - leakE) + Kg * powN * (Vn1 - Vn2 - Ke ) + NaG * powM * h * (Vn1 - Vn2 - NaE ));
250  return r;
251  }
252 
253  template <typename ScalarT>
254  static ScalarT kcl2EquQ( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& memC )
255  {
256  ScalarT r = -1.0 * memC * (Vn1 - Vn2);
257  return r;
258  }
259 
260  // n conservation equation
261  template <typename ScalarT>
262  static ScalarT nEquF( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& n, const ScalarT& Vrest )
263  {
264  ScalarT r = alphaN<ScalarT>( Vn1, Vn2, Vrest ) * (1.0 - n ) - betaN<ScalarT>( Vn1, Vn2, Vrest ) * n;
265  return r;
266  }
267 
268  template <typename ScalarT>
269  static ScalarT nEquQ( const ScalarT& n )
270  {
271  ScalarT r = -n;
272  return r;
273  }
274 
275  // m conservation equation
276  template <typename ScalarT>
277  static ScalarT mEquF( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& m, const ScalarT& Vrest )
278  {
279  ScalarT r = alphaM<ScalarT>( Vn1, Vn2, Vrest ) * (1.0 - m ) - betaM<ScalarT>( Vn1, Vn2, Vrest ) * m;
280  return r;
281  }
282 
283  template <typename ScalarT>
284  static ScalarT mEquQ( const ScalarT& m )
285  {
286  ScalarT r = -m;
287  return r;
288  }
289 
290  // h conservation equation
291  template <typename ScalarT>
292  static ScalarT hEquF( const ScalarT& Vn1, const ScalarT& Vn2, const ScalarT& h, const ScalarT& Vrest )
293  {
294  ScalarT r = alphaH<ScalarT>( Vn1, Vn2, Vrest ) * (1.0 - h ) - betaH<ScalarT>( Vn1, Vn2, Vrest ) * h;
295  return r;
296  }
297 
298  template <typename ScalarT>
299  static ScalarT hEquQ( const ScalarT& h )
300  {
301  ScalarT r = -h;
302  return r;
303  }
304 
305 public:
306  // iterator reference to the Neuron model which owns this instance.
307  // Getters and setters
309  {
310  return model_;
311  }
312 
313 private:
314 
315  Model & model_; //< Owning model
316 
317  // constant threshold adjustment
318  static const double VT;
319 
320  // derrived quantities computed in updateIntermediateVars
321  // and used in the load functions
332 
333  // state variables
336 
337  // local state indices (offsets)
340 
341  // local solution indices (offsets)
342  int li_Pos; // local index to positive node on this device
343  int li_Neg; // local index to negative node on this device
344  int li_nPro; // local index to n promoter value (Na current)
345  int li_mPro; // local index to m promoter value (K current)
346  int li_hPro; // local index to h promoter value (K current)
347 
348  // Matrix equation index variables:
349 
350  // Offset variables corresponding to the above declared indices.
356 
362 
366 
370 
374 };
375 
376 //-----------------------------------------------------------------------------
377 // Class : Model
378 // Purpose :
379 // Special Notes :
380 // Creator : Christy Warrender, SNL, Cognitive Modeling
381 // Creation Date : 06/22/12
382 //-----------------------------------------------------------------------------
383 class Model : public DeviceModel
384 {
385  typedef std::vector<Instance *> InstanceVector;
386 
387  friend class ParametricData<Model>;
388  friend class Instance;
389  friend class Traits;
390 
391 public:
392  Model(
393  const Configuration & configuration,
394  const ModelBlock & MB,
395  const FactoryBlock & factory_block);
396  ~Model();
397 
398 private:
399  Model();
400  Model(const Model &);
401  Model &operator=(const Model &);
402 
403 public:
404  virtual void forEachInstance(DeviceInstanceOp &op) const /* override */;
405 
406  virtual std::ostream &printOutInstances(std::ostream &os) const;
407 
408  bool processParams ();
409  bool processInstanceParams ();
410 
411 private:
412 
413  // parameter variables
414  double cMem; // membrane capacitance
415  double gMem; // membrane conductance of leak current
416  double eLeak; // reversal potential of leak current
417  double eNa; // sodium reversal potential
418  double gNa; // sodium base conductance
419  double eK; // potassium reversal potential
420  double gK; // potassium base conductance
421  double vRest; // resting potential
422 
423  // flags that parameters were given
424  bool cMemGiven;
425  bool gMemGiven;
427  bool eNaGiven;
428  bool gNaGiven;
429  bool eKGiven;
430  bool gKGiven;
432 
433 
434 public:
435  void addInstance(Instance *instance)
436  {
437  instanceContainer.push_back(instance);
438  }
439 
441  {
442  return instanceContainer;
443  }
444 
446  {
447  return instanceContainer;
448  }
449 
450 private:
451  std::vector<Instance*> instanceContainer;
452 };
453 
454 
455 //-----------------------------------------------------------------------------
456 // Class : Master
457 // Purpose :
458 // Special Notes :
459 // Creator : Christy Warrender, SNL, Cognitive Modeling
460 // Creation Date : 06/01/12
461 //-----------------------------------------------------------------------------
462 class Master : public DeviceMaster<Traits>
463 {
464 public:
466  const Configuration & configuration,
467  const FactoryBlock & factory_block,
468  const SolverState & ss1,
469  const DeviceOptions & do1)
470  : DeviceMaster<Traits>(configuration, factory_block, ss1, do1)
471  {}
472 
473  virtual bool updateState (double * solVec, double * staVec, double * stoVec);
474 };
475 
476 void registerDevice();
477 
478 } // namespace Neuron9
479 } // namespace Device
480 } // namespace Xyce
481 
485 
486 #endif