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