NEURON
impedanc.cpp
Go to the documentation of this file.
1 #include <../../nrnconf.h>
2 #undef check
3 #include "nrniv_mf.h"
4 #include "nrnmpi.h"
5 #include "nrn_ansi.h"
6 #include "nonlinz.h"
7 #include <InterViews/resource.h>
8 #include <complex>
9 #include "nrnoc2iv.h"
10 #include "classreg.h"
11 #include "membfunc.h"
12 
13 typedef void (*Pfrv4)(int, Node**, double**, Datum**);
14 
15 class Imp {
16  public:
17  Imp() = default;
18  virtual ~Imp();
19  // v(x)/i(x) and v(loc)/i(x) == v(x)/i(loc)
20  int compute(double freq, bool nonlin = false, int maxiter = 500);
21  void location(Section*, double);
22  double transfer_amp(Section*, double);
23  double input_amp(Section*, double);
24  double transfer_phase(Section*, double);
25  double input_phase(Section*, double);
26  double ratio_amp(Section*, double);
27 
28  private:
29  int loc(Section*, double);
30  void alloc();
31  void impfree();
32  void check();
33  void setmat(double);
34  void setmat1();
35 
36  void LUDecomp();
37  void solve();
38 
39  public:
40  double deltafac_ = .001;
41 
42  private:
43  int n = 0;
44  std::complex<double>* transfer = nullptr;
45  std::complex<double>* input = nullptr;
46  std::complex<double>* d = nullptr; /* diagonal */
47  std::complex<double>* pivot = nullptr;
48  int istim = -1; /* where current injected */
49  Section* sloc_ = nullptr;
50  double xloc_ = 0.;
51  NonLinImp* nli_ = nullptr;
52 };
53 
54 static void* cons(Object*) {
55  Imp* imp = new Imp();
56  return (void*) imp;
57 }
58 
59 static void destruct(void* v) {
60  Imp* imp = (Imp*) v;
61  delete imp;
62 }
63 
64 static double compute(void* v) {
65  Imp* imp = (Imp*) v;
66  int rval = 0;
67  bool nonlin = false;
68  if (ifarg(2)) {
69  nonlin = *getarg(2) ? true : false;
70  }
71  if (ifarg(3)) {
72  rval = imp->compute(*getarg(1), nonlin, int(chkarg(3, 1, 1e9)));
73  } else {
74  rval = imp->compute(*getarg(1), nonlin);
75  }
76  return double(rval);
77 }
78 
79 static double location(void* v) {
80  Imp* imp = (Imp*) v;
81  double x;
82  Section* sec = nullptr;
83  if (hoc_is_double_arg(1)) {
84  x = chkarg(1, -1., 1.);
85  if (x >= 0.0) {
86  sec = chk_access();
87  }
88  } else {
89  nrn_seg_or_x_arg(1, &sec, &x);
90  }
91  imp->location(sec, x);
92  return 0.;
93 }
94 
95 static double transfer_amp(void* v) {
96  Imp* imp = (Imp*) v;
97  double x;
98  Section* sec;
99  nrn_seg_or_x_arg(1, &sec, &x);
100  return imp->transfer_amp(sec, x);
101 }
102 static double input_amp(void* v) {
103  Imp* imp = (Imp*) v;
104  double x;
105  Section* sec;
106  nrn_seg_or_x_arg(1, &sec, &x);
107  return imp->input_amp(sec, x);
108 }
109 static double transfer_phase(void* v) {
110  Imp* imp = (Imp*) v;
111  double x;
112  Section* sec;
113  nrn_seg_or_x_arg(1, &sec, &x);
114  return imp->transfer_phase(sec, x);
115 }
116 static double input_phase(void* v) {
117  Imp* imp = (Imp*) v;
118  double x;
119  Section* sec;
120  nrn_seg_or_x_arg(1, &sec, &x);
121  return imp->input_phase(sec, x);
122 }
123 static double ratio_amp(void* v) {
124  Imp* imp = (Imp*) v;
125  double x;
126  Section* sec;
127  nrn_seg_or_x_arg(1, &sec, &x);
128  return imp->ratio_amp(sec, x);
129 }
130 
131 static double deltafac(void* v) {
132  Imp* imp = (Imp*) v;
133  if (ifarg(1)) {
134  imp->deltafac_ = chkarg(1, 1e-10, 1);
135  }
136  return imp->deltafac_;
137 }
138 
139 static Member_func members[] = {{"compute", compute},
140  {"loc", location},
141  {"input", input_amp},
142  {"transfer", transfer_amp},
143  {"ratio", ratio_amp},
144  {"input_phase", input_phase},
145  {"transfer_phase", transfer_phase},
146  {"deltafac", deltafac},
147  {nullptr, nullptr}};
148 
150  class2oc("Impedance", cons, destruct, members, nullptr, nullptr);
151 }
152 
154  if (sloc_) {
156  }
157  impfree();
158 }
159 
160 void Imp::impfree() {
161  if (d) {
162  delete[] d;
163  delete[] transfer;
164  delete[] input;
165  delete[] pivot;
166  d = nullptr;
167  }
168  if (nli_) {
169  delete nli_;
170  nli_ = nullptr;
171  }
172 }
173 
174 void Imp::check() {
175  NrnThread* _nt = nrn_threads;
176  nrn_thread_error("Impedance works with only one thread");
177  if (sloc_ && !sloc_->prop) {
179  sloc_ = nullptr;
180  }
181  if (tree_changed) {
182  setup_topology();
183  }
184  if (v_structure_change) {
185  recalc_diam();
186  }
187  if (n != _nt->end) {
188  alloc();
189  }
190 }
191 
192 void Imp::alloc() {
193  NrnThread* _nt = nrn_threads;
194  impfree();
195  n = _nt->end;
196  d = new std::complex<double>[n];
197  transfer = new std::complex<double>[n];
198  input = new std::complex<double>[n];
199  pivot = new std::complex<double>[n];
200 }
201 int Imp::loc(Section* sec, double x) {
202  if (x < 0.0 || sec == nullptr) {
203  return -1;
204  }
205  const Node* nd = node_exact(sec, x);
206  return nd->v_node_index;
207 }
208 
209 double Imp::transfer_amp(Section* sec, double x) {
210  check();
211  int vloc = loc(sec, x);
212  return nli_ ? nli_->transfer_amp(istim, vloc) : abs(transfer[vloc]);
213 }
214 
215 double Imp::input_amp(Section* sec, double x) {
216  check();
217  return nli_ ? nli_->input_amp(loc(sec, x)) : abs(input[loc(sec, x)]);
218 }
219 
220 double Imp::transfer_phase(Section* sec, double x) {
221  check();
222  return nli_ ? nli_->transfer_phase(istim, loc(sec, x)) : arg(transfer[loc(sec, x)]);
223 }
224 
225 double Imp::input_phase(Section* sec, double x) {
226  check();
227  return nli_ ? nli_->input_phase(loc(sec, x)) : arg(input[loc(sec, x)]);
228 }
229 
230 double Imp::ratio_amp(Section* sec, double x) {
231  check();
232  int i = loc(sec, x);
233  return nli_ ? nli_->ratio_amp(i, istim) : (abs(transfer[i] / input[i]));
234 }
235 
236 void Imp::location(Section* sec, double x) {
237  if (sloc_) {
239  }
240  sloc_ = sec;
241  xloc_ = x;
242  if (sloc_) {
244  }
245 }
246 
247 int Imp::compute(double freq, bool nonlin, int maxiter) {
248  int rval = 0;
249  check();
250  if (sloc_) {
251  istim = loc(sloc_, xloc_);
252  } else {
253  istim = -1;
254  if (nrnmpi_numprocs == 0) {
255  hoc_execerror("Impedance stimulus location is not specified.", 0);
256  }
257  }
258  if (n == 0 && nrnmpi_numprocs == 1)
259  return rval;
260  double omega = 1e-6 * 2 * 3.14159265358979323846 * freq; // wC has units of mho/cm2
261  if (nonlin) {
262  if (!nli_) {
263  nli_ = new NonLinImp();
264  }
265  nli_->compute(omega, deltafac_, maxiter);
266  rval = nli_->solve(istim);
267  } else {
268  if (nli_) {
269  delete nli_;
270  nli_ = nullptr;
271  }
272  if (istim == -1) {
273  hoc_execerror("Impedance stimulus location is not specified.", 0);
274  }
275  setmat(omega);
276  LUDecomp();
277  solve();
278  }
279  return rval;
280 }
281 
282 void Imp::setmat(double omega) {
283  const NrnThread* _nt = nrn_threads;
284  setmat1();
285  for (int i = 0; i < n; ++i) {
286  d[i] = std::complex<double>(NODED(_nt->_v_node[i]), NODERHS(_nt->_v_node[i]) * omega);
287  transfer[i] = 0.;
288  }
289  transfer[istim] = 1.e2 / NODEAREA(_nt->_v_node[istim]); // injecting 1nA
290  // rhs returned is then in units of mV or MegOhms
291 }
292 
293 
294 void Imp::setmat1() {
295  /*
296  The calculated g is good til someone else
297  changes something having to do with the matrix.
298  */
299  auto const sorted_token = nrn_ensure_model_data_are_sorted();
300  const NrnThread* _nt = nrn_threads;
301  const Memb_list* mlc = _nt->tml->ml;
302  assert(_nt->tml->index == CAP);
303  for (int i = 0; i < nrn_nthread; ++i) {
304  double cj = nrn_threads[i].cj;
305  nrn_threads[i].cj = 0;
306  // not useful except that many model description set g while computing i
307  nrn_rhs(sorted_token, nrn_threads[i]);
308  nrn_lhs(sorted_token, nrn_threads[i]);
309  nrn_threads[i].cj = cj;
310  }
311  for (int i = 0; i < n; ++i) {
312  NODERHS(_nt->_v_node[i]) = 0;
313  }
314  for (int i = 0; i < mlc->nodecount; ++i) {
315  NODERHS(mlc->nodelist[i]) = mlc->data(i, 0);
316  }
317 }
318 
320  const NrnThread* _nt = nrn_threads;
321  for (int i = _nt->end - 1; i >= _nt->ncell; --i) {
322  int ip = _nt->_v_parent[i]->v_node_index;
323  pivot[i] = NODEA(_nt->_v_node[i]) / d[i];
324  d[ip] -= pivot[i] * NODEB(_nt->_v_node[i]);
325  }
326 }
327 
328 void Imp::solve() {
329  for (int j = 0; j < nrn_nthread; ++j) {
330  const NrnThread* _nt = nrn_threads + j;
331  for (int i = istim; i >= _nt->ncell; --i) {
332  int ip = _nt->_v_parent[i]->v_node_index;
333  transfer[ip] -= transfer[i] * pivot[i];
334  }
335  for (int i = 0; i < _nt->ncell; ++i) {
336  transfer[i] /= d[i];
337  input[i] = 1. / d[i];
338  }
339  for (int i = _nt->ncell; i < _nt->end; ++i) {
340  int ip = _nt->_v_parent[i]->v_node_index;
341  transfer[i] -= NODEB(_nt->_v_node[i]) * transfer[ip];
342  transfer[i] /= d[i];
343  input[i] = (std::complex<double>(1) + input[ip] * pivot[i] * NODEB(_nt->_v_node[i])) /
344  d[i];
345  }
346  // take into account area
347  for (int i = _nt->ncell; i < _nt->end; ++i) {
348  input[i] *= 1e2 / NODEAREA(_nt->_v_node[i]);
349  }
350  }
351 }
Section * chk_access()
Definition: cabcode.cpp:449
void setup_topology(void)
Definition: cabcode.cpp:1635
int tree_changed
Definition: cabcode.cpp:51
Node * node_exact(Section *sec, double x)
like node_index but give proper node when x is 0 or 1 as well as in between
Definition: cabcode.cpp:1800
Definition: impedanc.cpp:15
Imp()=default
double input_phase(Section *, double)
Definition: impedanc.cpp:225
double xloc_
Definition: impedanc.cpp:50
double transfer_amp(Section *, double)
Definition: impedanc.cpp:209
std::complex< double > * d
Definition: impedanc.cpp:46
std::complex< double > * input
Definition: impedanc.cpp:45
int n
Definition: impedanc.cpp:43
double transfer_phase(Section *, double)
Definition: impedanc.cpp:220
std::complex< double > * transfer
Definition: impedanc.cpp:44
void alloc()
Definition: impedanc.cpp:192
void solve()
Definition: impedanc.cpp:328
int loc(Section *, double)
Definition: impedanc.cpp:201
virtual ~Imp()
Definition: impedanc.cpp:153
Section * sloc_
Definition: impedanc.cpp:49
std::complex< double > * pivot
Definition: impedanc.cpp:47
double deltafac_
Definition: impedanc.cpp:40
void setmat1()
Definition: impedanc.cpp:294
int istim
Definition: impedanc.cpp:48
void check()
Definition: impedanc.cpp:174
double ratio_amp(Section *, double)
Definition: impedanc.cpp:230
double input_amp(Section *, double)
Definition: impedanc.cpp:215
void impfree()
Definition: impedanc.cpp:160
void location(Section *, double)
Definition: impedanc.cpp:236
void LUDecomp()
Definition: impedanc.cpp:319
int compute(double freq, bool nonlin=false, int maxiter=500)
Definition: impedanc.cpp:247
void setmat(double)
Definition: impedanc.cpp:282
NonLinImp * nli_
Definition: impedanc.cpp:51
double transfer_amp(int curloc, int vloc)
Definition: nonlinz.cpp:61
int solve(int curloc)
Definition: nonlinz.cpp:182
void compute(double omega, double deltafac, int maxiter)
Definition: nonlinz.cpp:117
double ratio_amp(int clmploc, int vloc)
Definition: nonlinz.cpp:105
double input_phase(int curloc)
Definition: nonlinz.cpp:93
double input_amp(int curloc)
Definition: nonlinz.cpp:71
double transfer_phase(int curloc, int vloc)
Definition: nonlinz.cpp:83
void class2oc(const char *, ctor_f *cons, dtor_f *destruct, Member_func *, Member_ret_obj_func *, Member_ret_str_func *)
Definition: hoc_oop.cpp:1631
#define v
Definition: md1redef.h:11
#define sec
Definition: md1redef.h:20
#define i
Definition: md1redef.h:19
double chkarg(int, double low, double high)
Definition: code2.cpp:626
int hoc_is_double_arg(int narg)
Definition: code.cpp:864
#define assert(ex)
Definition: hocassrt.h:24
#define getarg
Definition: hocdec.h:17
static double transfer_amp(void *v)
Definition: impedanc.cpp:95
void Impedance_reg()
Definition: impedanc.cpp:149
static double location(void *v)
Definition: impedanc.cpp:79
void(* Pfrv4)(int, Node **, double **, Datum **)
Definition: impedanc.cpp:13
static double input_amp(void *v)
Definition: impedanc.cpp:102
static Member_func members[]
Definition: impedanc.cpp:139
static double deltafac(void *v)
Definition: impedanc.cpp:131
static double input_phase(void *v)
Definition: impedanc.cpp:116
static void * cons(Object *)
Definition: impedanc.cpp:54
static void destruct(void *v)
Definition: impedanc.cpp:59
static double ratio_amp(void *v)
Definition: impedanc.cpp:123
static double compute(void *v)
Definition: impedanc.cpp:64
static double transfer_phase(void *v)
Definition: impedanc.cpp:109
#define CAP
Definition: membfunc.hpp:60
NrnThread * nrn_threads
Definition: multicore.cpp:56
int nrn_nthread
Definition: multicore.cpp:55
static void nrn_lhs(NrnThread *_nt)
int v_structure_change
Definition: nrnoc_aux.cpp:20
void hoc_execerror(const char *s1, const char *s2)
Definition: nrnoc_aux.cpp:39
static void nrn_rhs(NrnThread *_nt)
void section_ref(Section *)
Definition: solve.cpp:520
void section_unref(Section *)
Definition: solve.cpp:509
void nrn_seg_or_x_arg(int iarg, Section **psec, double *px)
Definition: point.cpp:170
void recalc_diam(void)
Definition: treeset.cpp:923
neuron::model_sorted_token nrn_ensure_model_data_are_sorted()
Ensure neuron::container::* data are sorted.
Definition: treeset.cpp:2182
size_t j
int ifarg(int)
Definition: code.cpp:1607
void nrn_thread_error(const char *s)
Definition: multicore.cpp:297
#define NODEAREA(n)
Definition: section.h:101
#define NODEA(n)
Definition: section_fwd.hpp:52
#define NODEB(n)
Definition: section_fwd.hpp:53
#define NODERHS(n)
Definition: section_fwd.hpp:55
#define NODED(n)
Definition: section_fwd.hpp:54
A view into a set of mechanism instances.
Definition: nrnoc_ml.h:34
int nodecount
Definition: nrnoc_ml.h:78
Node ** nodelist
Definition: nrnoc_ml.h:68
std::vector< double * > data()
Get a vector of double* representing the model data.
Definition: memblist.cpp:64
Definition: section.h:105
int v_node_index
Definition: section.h:212
Represent main neuron object computed by single thread.
Definition: multicore.h:58
NrnThreadMembList * tml
Definition: multicore.h:62
int ncell
Definition: multicore.h:64
int end
Definition: multicore.h:65
Node ** _v_parent
Definition: multicore.h:91
Node ** _v_node
Definition: multicore.h:90
Memb_list * ml
Definition: multicore.h:35
Definition: hocdec.h:173
Prop * prop
Definition: section.h:71
Non-template stable handle to a generic value.