Qucs-core  0.0.18
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
hbsolver.cpp
Go to the documentation of this file.
1 /*
2  * hbsolver.cpp - harmonic balance solver class implementation
3  *
4  * Copyright (C) 2005, 2006, 2007, 2008 Stefan Jahn <stefan@lkcc.org>
5  *
6  * This is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2, or (at your option)
9  * any later version.
10  *
11  * This software is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this package; see the file COPYING. If not, write to
18  * the Free Software Foundation, Inc., 51 Franklin Street - Fifth Floor,
19  * Boston, MA 02110-1301, USA.
20  *
21  * $Id$
22  *
23  */
24 
25 #if HAVE_CONFIG_H
26 # include <config.h>
27 #endif
28 
29 #include <stdio.h>
30 
31 #include "object.h"
32 #include "logging.h"
33 #include "complex.h"
34 #include "vector.h"
35 #include "node.h"
36 #include "circuit.h"
37 #include "component_id.h"
38 #include "net.h"
39 #include "netdefs.h"
40 #include "strlist.h"
41 #include "ptrlist.h"
42 #include "tvector.h"
43 #include "tmatrix.h"
44 #include "eqnsys.h"
45 #include "analysis.h"
46 #include "dataset.h"
47 #include "fourier.h"
48 #include "hbsolver.h"
49 
50 #define HB_DEBUG 0
51 
52 namespace qucs {
53 
54 using namespace fourier;
55 
56 // Constructor creates an unnamed instance of the hbsolver class.
57 hbsolver::hbsolver () : analysis () {
59  frequency = 0;
60  nlnodes = lnnodes = banodes = nanodes = NULL;
61  Y = Z = A = NULL;
62  NA = YV = JQ = JG = JF = NULL;
63  OM = IR = QR = RH = IG = FQ = VS = VP = FV = IL = IN = IC = IS = NULL;
64  vs = x = NULL;
65  runs = 0;
66  ndfreqs = NULL;
67 }
68 
69 // Constructor creates a named instance of the hbsolver class.
70 hbsolver::hbsolver (char * n) : analysis (n) {
72  frequency = 0;
73  nlnodes = lnnodes = banodes = nanodes = NULL;
74  Y = Z = A = NULL;
75  NA = YV = JQ = JG = JF = NULL;
76  OM = IR = QR = RH = IG = FQ = VS = VP = FV = IL = IN = IC = IS = NULL;
77  vs = x = NULL;
78  runs = 0;
79  ndfreqs = NULL;
80 }
81 
82 // Destructor deletes the hbsolver class object.
84  // delete nodelists
85  if (nlnodes) delete nlnodes;
86  if (lnnodes) delete lnnodes;
87  if (banodes) delete banodes;
88  if (nanodes) delete nanodes;
89 
90  // delete temporary matrices
91  if (A) delete A;
92  if (Z) delete Z;
93  if (Y) delete Y;
94 
95  // delete matrices
96  if (NA) delete NA;
97  if (YV) delete YV;
98  if (JQ) delete JQ;
99  if (JG) delete JG;
100  if (JF) delete JF;
101 
102  // delete vectors
103  if (IC) delete IC;
104  if (IS) delete IS;
105  if (FV) delete FV;
106  if (IL) delete IL;
107  if (IN) delete IN;
108  if (IG) delete IG;
109  if (FQ) delete FQ;
110  if (VS) delete VS;
111  if (VP) delete VP;
112  if (vs) delete vs;
113  if (OM) delete OM;
114  if (IR) delete IR;
115  if (QR) delete QR;
116  if (RH) delete RH;
117 
118  if (x) delete x;
119  if (ndfreqs) delete[] ndfreqs;
120 }
121 
122 /* The copy constructor creates a new instance of the hbsolver class
123  based on the given hbsolver object. */
124 hbsolver::hbsolver (hbsolver & o) : analysis (o) {
125  frequency = o.frequency;
126  negfreqs = o.negfreqs;
127  posfreqs = o.posfreqs;
128  nlnodes = o.nlnodes;
129  lnnodes = o.lnnodes;
130  banodes = o.banodes;
131  nanodes = o.nanodes;
132  Y = Z = A = NULL;
133  NA = YV = JQ = JG = JF = NULL;
134  OM = IR = QR = RH = IG = FQ = VS = VP = FV = IL = IN = IC = IS = NULL;
135  vs = x = NULL;
136  runs = o.runs;
137  ndfreqs = NULL;
138 }
139 
140 #define VS_(r) (*VS) (r)
141 #define OM_(r) (*OM) (r)
142 
143 /* This is the HB netlist solver. It prepares the circuit list and
144  solves it then. */
145 int hbsolver::solve (void) {
146 
147  int iterations = 0, done = 0;
148  int MaxIterations = getPropertyInteger ("MaxIter");
149 
150  // collect different parts of the circuit
151  splitCircuits ();
152 
153  // create frequency array
155 
156  // find interconnects between the linear and non-linear subcircuit
157  getNodeLists ();
158 
159  // prepares the linear part --> 0 = IC + [YV] * VS
160  prepareLinear ();
161 
162  runs++;
163  logprint (LOG_STATUS, "NOTIFY: %s: solving for %d frequencies\n",
164  getName (), lnfreqs);
165 
166  if (nbanodes > 0) {
167 
168  // start balancing
169  logprint (LOG_STATUS, "NOTIFY: %s: balancing at %d nodes\n", getName (),
170  nbanodes);
171 
172  // prepares the non-linear part
173  prepareNonLinear ();
174 
175 #if HB_DEBUG
176  fprintf (stderr, "YV -- transY in f:\n"); YV->print ();
177  fprintf (stderr, "IC -- constant current in f:\n"); IC->print ();
178 #endif
179 
180  // start iteration
181  do {
182  iterations++;
183 
184 #if HB_DEBUG
185  fprintf (stderr, "\n -- iteration step: %d\n", iterations);
186  fprintf (stderr, "vs -- voltage in t:\n"); vs->print ();
187 #endif
188 
189  // evaluate component functionality and fill matrices and vectors
190  loadMatrices ();
191 
192 #if HB_DEBUG
193  fprintf (stderr, "FQ -- charge in t:\n"); FQ->print ();
194  fprintf (stderr, "IG -- current in t:\n"); IG->print ();
195 #endif
196 
197  // currents into frequency domain
198  VectorFFT (IG);
199 
200  // charges into frequency domain
201  VectorFFT (FQ);
202 
203  // right hand side currents and charges into the frequency domain
204  VectorFFT (IR);
205  VectorFFT (QR);
206 
207 #if HB_DEBUG
208  fprintf (stderr, "VS -- voltage in f:\n"); VS->print ();
209  fprintf (stderr, "FQ -- charge in f:\n"); FQ->print ();
210  fprintf (stderr, "IG -- current in f:\n"); IG->print ();
211  fprintf (stderr, "IR -- corrected Jacobi current in f:\n"); IR->print ();
212 #endif
213 
214  // solve HB equation --> FV = IC + [YV] * VS + j[O] * FQ + IG
215  solveHB ();
216 
217 #if HB_DEBUG
218  fprintf (stderr, "FV -- error vector F(V) in f:\n"); FV->print ();
219  fprintf (stderr, "IL -- linear currents in f:\n"); IL->print ();
220  fprintf (stderr, "IN -- non-linear currents in f:\n"); IN->print ();
221  fprintf (stderr, "RH -- right-hand side currents in f:\n"); RH->print ();
222 #endif
223 
224  // termination criteria met
225  if (iterations > 1 && checkBalance ()) {
226  done = 1;
227  break;
228  }
229 
230 #if HB_DEBUG
231  fprintf (stderr, "JG -- G-Jacobian in t:\n"); JG->print ();
232  fprintf (stderr, "JQ -- C-Jacobian in t:\n"); JQ->print ();
233 #endif
234 
235  // G-Jacobian into frequency domain
236  MatrixFFT (JG);
237 
238  // C-Jacobian into frequency domain
239  MatrixFFT (JQ);
240 
241 #if HB_DEBUG
242  fprintf (stderr, "JQ -- dQ/dV C-Jacobian in f:\n"); JQ->print ();
243  fprintf (stderr, "JG -- dI/dV G-Jacobian in f:\n"); JG->print ();
244 #endif
245 
246  // calculate Jacobian --> JF = [YV] + j[O] * JQ + JG
247  calcJacobian ();
248 
249 #if HB_DEBUG
250  fprintf (stderr, "JF -- full Jacobian in f:\n"); JF->print ();
251 #endif
252 
253  // solve equation system --> JF * VS(n+1) = JF * VS(n) - FV
254  solveVoltages ();
255 
256 #if HB_DEBUG
257  fprintf (stderr, "VS -- next voltage in f:\n"); VS->print ();
258 #endif
259 
260  // inverse FFT of frequency domain voltage vector VS(n+1)
261  VectorIFFT (vs);
262  }
263  // check termination criteria (balanced frequency domain currents)
264  while (!done && iterations < MaxIterations);
265 
266  if (iterations >= MaxIterations) {
268  e->setText ("no convergence in %s analysis after %d iterations",
269  getName (), iterations);
270  throw_exception (e);
271  logprint (LOG_ERROR, "%s: no convergence after %d iterations\n",
272  getName (), iterations);
273  }
274  else {
275  logprint (LOG_STATUS, "%s: convergence reached after %d iterations\n",
276  getName (), iterations);
277  }
278  }
279  else {
280  // no balancing necessary
281  logprint (LOG_STATUS, "NOTIFY: %s: no balancing necessary\n", getName ());
282  }
283 
284  // print exception stack
285  estack.print ();
286 
287  // apply AC analysis to the complete network in order to obtain the
288  // final results
289  finalSolution ();
290 
291  // save results into dataset
292  saveResults ();
293  return 0;
294 }
295 
296 /* Goes through the list of circuit objects and runs its calcHB()
297  function. */
298 void hbsolver::calc (hbsolver * self) {
299  circuit * root = self->getNet()->getRoot ();
300  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ()) {
301  c->calcHB (self->frequency);
302  }
303 }
304 
305 /* Goes through the list of circuit objects and runs its initHB()
306  function. */
307 void hbsolver::initHB (void) {
308  circuit * root = subnet->getRoot ();
309  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ()) {
310  c->initHB ();
311  }
312 }
313 
314 /* Goes through the list of circuit objects and runs its initDC()
315  function. */
316 void hbsolver::initDC (void) {
317  circuit * root = subnet->getRoot ();
318  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ()) {
319  c->initDC ();
320  }
321 }
322 
323 // Returns true if circuit is a HB source.
325  int type = c->getType ();
326  if (type == CIR_PAC || type == CIR_VAC || type == CIR_VDC)
327  return true;
328  return false;
329 }
330 
331 // Expands the frequency array using the given frequency and the order.
332 void hbsolver::expandFrequencies (nr_double_t f, int n) {
335  int i, k, len = nfreqs.getSize ();
336  negfreqs.clear ();
337  posfreqs.clear ();
338  if (len > 0) {
339  // frequency expansion for full frequency sets
340  for (i = 0; i <= n + 1; i++) {
341  for (k = 0; k < len; k++) {
342  negfreqs.add (i * f + nfreqs.get (k));
343  }
344  }
345  for (i = -n; i < 0; i++) {
346  for (k = 0; k < len; k++) {
347  negfreqs.add (i * f + nfreqs.get (k));
348  }
349  }
350  for (i = 0; i <= 2 * n + 1; i++) {
351  for (k = 0; k < len; k++) {
352  posfreqs.add (i * f + pfreqs.get (k));
353  }
354  }
355  }
356  else {
357  // first frequency
358  for (i = 0; i <= n + 1; i++) negfreqs.add (i * f);
359  for (i = -n; i < 0; i++) negfreqs.add (i * f);
360  for (i = 0; i <= 2 * n + 1; i++) posfreqs.add (i * f);
361  }
362 }
363 
364 // Calculates an order fulfilling the "power of two" requirement.
366  int o, order = n * 2; // current order + DC + negative freqencies
367  for (o = 1; o < order; o <<= 1) ; // a power of 2
368  return o / 2 - 1;
369 }
370 
371 /* The function computes the harmonic frequencies excited in the
372  circuit list depending on the maximum number of harmonics per
373  exitation and saves its results into the 'negfreqs' vector. */
375 
376  // initialization
377  negfreqs.clear ();
378  posfreqs.clear ();
379  rfreqs.clear ();
380  dfreqs.clear ();
381  if (ndfreqs) delete[] ndfreqs;
382 
383  // obtain order
384  int i, n = calcOrder (getPropertyInteger ("n"));
385 
386  // expand frequencies for each exitation
387  nr_double_t f;
388  for (ptrlistiterator<circuit> it (excitations); *it; ++it) {
389  circuit * c = it.current ();
390  if (c->getType () != CIR_VDC) { // no extra DC sources
391  if ((f = c->getPropertyDouble ("f")) != 0.0) {
392  if (!dfreqs.contains (f)) { // no double frequencies
393  dfreqs.add (f);
394  expandFrequencies (f, n);
395  }
396  }
397  }
398  }
399 
400  // no excitations
401  if (negfreqs.getSize () == 0) {
402  // use specified frequency
403  f = getPropertyDouble ("f");
404  dfreqs.add (f);
405  expandFrequencies (f, n);
406  }
407 
408  // build frequency dimension lengths
409  ndfreqs = new int[dfreqs.getSize ()];
410  for (i = 0; i < dfreqs.getSize (); i++) {
411  ndfreqs[i] = (n + 1) * 2;
412  }
413 
414 #if HB_DEBUG
415  fprintf (stderr, "%d frequencies: [ ", negfreqs.getSize ());
416  for (i = 0; i < negfreqs.getSize (); i++) {
417  fprintf (stderr, "%g ", (double) real (negfreqs.get (i)));
418  }
419  fprintf (stderr, "]\n");
420 #endif /* HB_DEBUG */
421 
422  // build list of positive frequencies including DC
423  for (n = 0; n < negfreqs.getSize (); n++) {
424  if ((f = negfreqs (n)) < 0.0) continue;
425  rfreqs.add (f);
426  }
427  lnfreqs = rfreqs.getSize ();
428  nlfreqs = negfreqs.getSize ();
429 
430  // pre-calculate the j[O] vector
432  for (n = i = 0; n < nlfreqs; n++, i++)
433  OM_(n) = nr_complex_t (0, 2 * M_PI * negfreqs (i));
434 }
435 
436 // Split netlist into excitation, linear and non-linear part.
438  circuit * root = subnet->getRoot ();
439  for (circuit * c = root; c != NULL; c = (circuit *) c->getNext ()) {
440  if (c->isNonLinear ()) {
441  // non-linear part
442  nolcircuits.add (c);
443  }
444  else if (isExcitation (c)) {
445  // get sinusoidal sources
446  excitations.add (c);
447  }
448  else if (c->getType () != CIR_GROUND) {
449  // linear part
450  lincircuits.add (c);
451  }
452  }
453 }
454 
455 // Creates a nodelist for the given circuit list.
456 strlist * hbsolver::circuitNodes (ptrlist<circuit> circuits) {
457  strlist * nodes = new strlist ();
458  for (ptrlistiterator<circuit> it (circuits); *it; ++it) {
459  circuit * c = it.current ();
460  for (int i = 0; i < c->getSize (); i++) {
461  char * n = c->getNode(i)->getName ();
462  if (strcmp (n, "gnd")) {
463  if (!nodes->contains (n)) nodes->add (n);
464  }
465  }
466  }
467  return nodes;
468 }
469 
470 // Obtain node lists for linear and non-linear part.
472  // non-linear nodes
474  // linear nodes
476  // excitation nodes
478 
479 #if DEBUG && 0
480  fprintf (stderr, "nonlinear nodes: [ %s ]\n", nlnodes->toString ());
481  fprintf (stderr, " linear nodes: [ %s ]\n", lnnodes->toString ());
482 #endif
483 
484  // organization of the nodes for the MNA:
485  // --------------------------------------
486  // 1.) balanced nodes: all connected to at least one non-linear device
487  // 2.a) the excitation nodes
488  // 2.b) all linear nodes not contained in 1. and 2.a
489  // 2.c) additional gyrators of linear nodes (built-in voltage sources)
490  // please note: excitation nodes also in 2.b; 1. and 2.a are 'ports'
491 
492  nanodes = new strlist (*nlnodes); // list 1.
493  strlistiterator it;
494  // add excitation nodes; list 2.a
495  for (it = strlistiterator (exnodes); *it; ++it)
496  nanodes->append (*it);
497  // add linear nodes; list 2.b
498  for (it = strlistiterator (lnnodes); *it; ++it) {
499  if (!nanodes->contains (*it))
500  nanodes->append (*it);
501  }
502 
503  banodes = new strlist (*nlnodes);
504 #if DEBUG && 0
505  fprintf (stderr, " balanced nodes: [ %s ]\n", banodes->toString ());
506  fprintf (stderr, " exnodes nodes: [ %s ]\n", exnodes->toString ());
507  fprintf (stderr, " nanodes nodes: [ %s ]\n", nanodes->toString ());
508 #endif
509 }
510 
511 /* The function applies unique voltage source identifiers to each
512  built in internal voltage source in the list of circuits. */
513 int hbsolver::assignVoltageSources (ptrlist<circuit> circuits) {
514  int sources = 0;
515  for (ptrlistiterator<circuit> it (circuits); *it; ++it) {
516  circuit * c = it.current ();
517  if (c->getVoltageSources () > 0) {
518  c->setVoltageSource (sources);
519  sources += c->getVoltageSources ();
520  }
521  }
522  return sources;
523 }
524 
525 /* The function assigns unique node identifiers to each circuit node. */
526 int hbsolver::assignNodes (ptrlist<circuit> circuits, strlist * nodes,
527  int offset) {
528  // through all collected nodes
529  for (int nr = 0; nr < nodes->length (); nr++) {
530  char * nn = nodes->get (nr);
531  // through all circuits
532  for (ptrlistiterator<circuit> it (circuits); *it; ++it) {
533  circuit * c = it.current ();
534  // assign current identifier if any of the circuit nodes equals
535  // the current one
536  for (int i = 0; i < c->getSize (); i++) {
537  node * n = c->getNode (i);
538  if (!strcmp (n->getName (), nn)) n->setNode (offset + nr + 1);
539  }
540  }
541  }
542  return nodes->length ();
543 }
544 
545 // Prepares the linear operations.
547  for (ptrlistiterator<circuit> it (lincircuits); *it; ++it) (*it)->initHB ();
549  nnlvsrcs = excitations.length ();
550  nnanodes = nanodes->length ();
551  nexnodes = exnodes->length ();
552  nbanodes = banodes->length ();
558 }
559 
560 /* The function creates the complex linear network MNA matrix. It
561  contains the MNA entries for all linear components for each
562  requested frequency. */
564  int M = nlnvsrcs;
565  int N = nnanodes;
566  int f = 0;
567  nr_double_t freq;
568 
569  // create new MNA matrix
570  A = new tmatrix<nr_complex_t> ((N + M) * lnfreqs);
571 
572  // through each frequency
573  for (int i = 0; i < rfreqs.getSize (); i++) {
574  freq = rfreqs (i);
575  // calculate components' MNA matrix for the given frequency
576  for (ptrlistiterator<circuit> it (lincircuits); *it; ++it)
577  (*it)->calcHB (freq);
578  // fill in all matrix entries for the given frequency
579  fillMatrixLinearA (A, f++);
580  }
581 
582  // save a copy of the original MNA matrix
583  NA = new tmatrix<nr_complex_t> (*A);
584 }
585 
586 // some definitions for the linear matrix filler
587 #undef A_
588 #undef B_
589 #define A_(r,c) (*A) ((r)*lnfreqs+f,(c)*lnfreqs+f)
590 #define G_(r,c) A_(r,c)
591 #define B_(r,c) A_(r,c+N)
592 #define C_(r,c) A_(r+N,c)
593 #define D_(r,c) A_(r+N,c+N)
594 
595 /* This function fills in the MNA matrix entries into the A matrix for
596  a given frequency index. */
598  int N = nnanodes;
599 
600  // through each linear circuit
601  for (ptrlistiterator<circuit> it (lincircuits); *it; ++it) {
602  circuit * cir = it.current ();
603  int s = cir->getSize ();
604  int nr, nc, r, c, v;
605 
606  // apply G-matrix entries
607  for (r = 0; r < s; r++) {
608  if ((nr = cir->getNode(r)->getNode () - 1) < 0) continue;
609  for (c = 0; c < s; c++) {
610  if ((nc = cir->getNode(c)->getNode () - 1) < 0) continue;
611  G_(nr, nc) += cir->getY (r, c);
612  }
613  }
614 
615  // augmented part -- built in voltage sources
616  if ((v = cir->getVoltageSources ()) > 0) {
617 
618  // apply B-matrix entries
619  for (r = 0; r < s; r++) {
620  if ((nr = cir->getNode(r)->getNode () - 1) < 0) continue;
621  for (c = 0; c < v; c++) {
622  nc = cir->getVoltageSource () + c;
623  B_(nr, nc) += cir->getB (r, nc);
624  }
625  }
626 
627  // apply C-matrix entries
628  for (r = 0; r < v; r++) {
629  nr = cir->getVoltageSource () + r;
630  for (c = 0; c < s; c++) {
631  if ((nc = cir->getNode(c)->getNode () - 1) < 0) continue;
632  C_(nr, nc) += cir->getC (nr, c);
633  }
634  }
635 
636  // apply D-matrix entries
637  for (r = 0; r < v; r++) {
638  nr = cir->getVoltageSource () + r;
639  for (c = 0; c < v; c++) {
640  nc = cir->getVoltageSource () + c;
641  D_(nr, nc) += cir->getD (nr, nc);
642  }
643  }
644  }
645  }
646 }
647 
648 // The function inverts the given matrix A into the matrix H.
650  tmatrix<nr_complex_t> * H) {
652  int N = A->getCols ();
655 
656  try_running () {
657  // create LU decomposition of the A matrix
659  eqns.passEquationSys (A, x, z);
660  eqns.solve ();
661  }
662  // appropriate exception handling
663  catch_exception () {
664  case EXCEPTION_PIVOT:
665  default:
666  logprint (LOG_ERROR, "WARNING: %s: during TI inversion\n", getName ());
667  estack.print ();
668  }
669 
670  // use the LU decomposition to obtain the inverse H
672  for (int c = 0; c < N; c++) {
673  z->set (0.0);
674  z->set (c, 1.0);
675  eqns.passEquationSys (A, x, z);
676  eqns.solve ();
677  for (int r = 0; r < N; r++) H->set (r, c, x->get (r));
678  }
679  delete x;
680  delete z;
681 }
682 
683 // Some defines for matrix element access.
684 #define V_(r) (*V) (r)
685 #define I_(r) (*I) (r)
686 #undef A_
687 #define A_(r,c) (*A) (r,c)
688 
689 #define Z_(r,c) (*Z) (r,c)
690 #define Y_(r,c) (*Y) (r,c)
691 
692 #define ZVU_(r,c) Z_(r,c)
693 #define ZVL_(r,c) Z_((r)*lnfreqs+f+sn,c)
694 #define ZCU_(r,c) Z_(r,(c)*lnfreqs+f+sn)
695 #define ZCL_(r,c) Z_((r)*lnfreqs+f+sn,(c)*lnfreqs+f+sn)
696 
697 #define YV_(r,c) (*YV) (r,c)
698 #define NA_(r,c) (*NA) (r,c)
699 #define JF_(r,c) (*JF) (r,c)
700 
701 /* The following function performs the following steps:
702  1. form the MNA matrix A including all nodes (linear, non-linear and
703  excitations)
704  2. compute the variable transimpedance matrix entries for the nodes
705  to be balanced
706  3. compute the constant transimpedance matrix entries for the constant
707  current vector caused by the excitations
708  4. invert this overall transimpedance matrix
709  5. extract the variable transadmittance matrix entries
710 */
712  int M = nlnvsrcs;
713  int N = nnanodes;
714  int c, r, f;
716 
717  // size of overall MNA matrix
718  int sa = (N + M) * lnfreqs;
719  int sv = nbanodes;
720  int se = nnlvsrcs;
721  int sy = sv + se;
722 
723  // allocate new transimpedance matrix
724  Z = new tmatrix<nr_complex_t> (sy * lnfreqs);
725 
726  // prepare equation system
730 
731  // 1. create variable transimpedance matrix entries relating
732  // voltages at the balanced nodes to the currents through these
733  // nodes into the non-linear part
734  int sn = sv * lnfreqs;
735  V = new tvector<nr_complex_t> (sa);
736  I = new tvector<nr_complex_t> (sa);
737 
738  // connect a 100 Ohm resistor (to ground) to balanced node in the MNA matrix
739  for (c = 0; c < sv * lnfreqs; c++) A_(c, c) += 0.01;
740 
741  // connect a 100 Ohm resistor (in parallel) to each excitation
742  for (ite = ptrlistiterator<circuit> (excitations); *ite; ++ite) {
743  circuit * vs = ite.current ();
744  // get positive and negative node
745  int pnode = vs->getNode(NODE_1)->getNode ();
746  int nnode = vs->getNode(NODE_2)->getNode ();
747  for (f = 0; f < lnfreqs; f++) { // for each frequency
748  int pn = (pnode - 1) * lnfreqs + f;
749  int nn = (nnode - 1) * lnfreqs + f;
750  if (pnode) A_(pn, pn) += 0.01;
751  if (nnode) A_(nn, nn) += 0.01;
752  if (pnode && nnode) {
753  A_(pn, nn) -= 0.01;
754  A_(nn, pn) -= 0.01;
755  }
756  }
757  }
758 
759  // LU decompose the MNA matrix
760  try_running () {
761  eqns.setAlgo (ALGO_LU_FACTORIZATION_CROUT);
762  eqns.passEquationSys (A, V, I);
763  eqns.solve ();
764  }
765  // appropriate exception handling
766  catch_exception () {
767  case EXCEPTION_PIVOT:
768  default:
769  logprint (LOG_ERROR, "WARNING: %s: during A factorization\n", getName ());
770  estack.print ();
771  }
772 
773  // aquire variable transimpedance matrix entries
774  eqns.setAlgo (ALGO_LU_SUBSTITUTION_CROUT);
775  for (c = 0; c < sn; c++) {
776  I->set (0.0);
777  I_(c) = 1.0;
778  eqns.passEquationSys (A, V, I);
779  eqns.solve ();
780  // ZV | ..
781  // ---+---
782  // .. | ..
783  for (r = 0; r < sn; r++) ZVU_(r, c) = V_(r);
784  // .. | ..
785  // ---+---
786  // ZV | ..
787  r = 0;
788  for (ite = ptrlistiterator<circuit> (excitations); *ite; ++ite, r++) {
789  // lower part entries
790  for (f = 0; f < lnfreqs; f++) {
791  ZVL_(r, c) = excitationZ (V, *ite, f);
792  }
793  }
794  }
795 
796  // create constant transimpedance matrix entries relating the
797  // source voltages to the interconnection currents
798  int vsrc = 0;
799  // aquire constant transadmittance matrix entries
800  for (ptrlistiterator<circuit> it (excitations); *it; ++it, vsrc++) {
801  circuit * vs = it.current ();
802  // get positive and negative node
803  int pnode = vs->getNode(NODE_1)->getNode ();
804  int nnode = vs->getNode(NODE_2)->getNode ();
805  for (f = 0; f < lnfreqs; f++) { // for each frequency
806  int pn = (pnode - 1) * lnfreqs + f;
807  int nn = (nnode - 1) * lnfreqs + f;
808  I->set (0.0);
809  if (pnode) I_(pn) = +1.0;
810  if (nnode) I_(nn) = -1.0;
811  eqns.passEquationSys (A, V, I);
812  eqns.solve ();
813  // .. | ZC
814  // ---+---
815  // .. | ..
816  for (r = 0; r < sn; r++) {
817  // upper part of the entries
818  ZCU_(r, vsrc) = V_(r);
819  }
820  // .. | ..
821  // ---+---
822  // .. | ZC
823  r = 0;
824  for (ite = ptrlistiterator<circuit> (excitations); *ite; ++ite, r++) {
825  // lower part entries
826  ZCL_(r, vsrc) = excitationZ (V, *ite, f);
827  }
828  }
829  }
830  delete I;
831  delete V;
832 
833  // allocate new transadmittance matrix
834  Y = new tmatrix<nr_complex_t> (sy * lnfreqs);
835 
836  // invert the Z matrix to a Y matrix
837  invertMatrix (Z, Y);
838 
839  // substract the 100 Ohm resistor
840  for (c = 0; c < sy * lnfreqs; c++) Y_(c, c) -= 0.01;
841 
842  // extract the variable transadmittance matrix
843  YV = new tmatrix<nr_complex_t> (sv * nlfreqs);
844 
845  // variable transadmittance matrix must be continued conjugately
846  *YV = expandMatrix (*Y, sv);
847 
848  // delete overall temporary MNA matrix
849  delete A; A = NULL;
850  // delete transimpedance matrix
851  delete Z; Z = NULL;
852 }
853 
854 /* Little helper function obtaining a transimpedance value for the
855  given voltage source (excitation) and for a given frequency
856  index. */
858  int f) {
859  // get positive and negative node
860  int pnode = vs->getNode(NODE_1)->getNode ();
861  int nnode = vs->getNode(NODE_2)->getNode ();
862  nr_complex_t z = 0.0;
863  if (pnode) z += V_((pnode - 1) * lnfreqs + f);
864  if (nnode) z -= V_((nnode - 1) * lnfreqs + f);
865  return z;
866 }
867 
868 /* This function computes the constant current vectors using the
869  voltage of the excitations and the transadmittance matrix
870  entries. */
872  int se = nnlvsrcs * lnfreqs;
873  int sn = nbanodes * lnfreqs;
874  int r, c, vsrc = 0;
875 
876  // collect excitation voltages
877  tvector<nr_complex_t> VC (se);
878  for (ptrlistiterator<circuit> it (excitations); *it; ++it, vsrc++) {
879  circuit * vs = it.current ();
880  vs->initHB ();
881  vs->setVoltageSource (0);
882  for (int f = 0; f < rfreqs.getSize (); f++) { // for each frequency
883  nr_double_t freq = rfreqs (f);
884  vs->calcHB (freq);
885  VC (vsrc * lnfreqs + f) = vs->getE (VSRC_1);
886  }
887  }
888 
889  // compute constant current vector for balanced nodes
890  IC = new tvector<nr_complex_t> (sn);
891  // .. | YC * VC
892  // ---+---
893  // .. | ..
894  for (r = 0; r < sn; r++) {
895  nr_complex_t i = 0.0;
896  for (c = 0; c < se; c++) {
897  i += Y_(r, c + sn) * VC (c);
898  }
899  int f = r % lnfreqs;
900  if (f != 0 && f != lnfreqs - 1) i /= 2;
901  IC->set (r, i);
902  }
903  // expand the constant current conjugate
904  *IC = expandVector (*IC, nbanodes);
905 
906  // compute constant current vector for sources itself
907  IS = new tvector<nr_complex_t> (se);
908  // .. | ..
909  // ---+---
910  // .. | YC * VC
911  for (r = 0; r < se; r++) {
912  nr_complex_t i = 0.0;
913  for (c = 0; c < se; c++) {
914  i += Y_(r + sn, c + sn) * VC (c);
915  }
916  IS->set (r, i);
917  }
918 
919  // delete overall transadmittance matrix
920  delete Y; Y = NULL;
921 }
922 
923 /* Checks whether currents through the interconnects of the linear and
924  non-linear subcircuit (in the frequency domain) are equal. */
926  nr_double_t iabstol = getPropertyDouble ("iabstol");
927  nr_double_t vabstol = getPropertyDouble ("vabstol");
928  nr_double_t reltol = getPropertyDouble ("reltol");
929  int n, len = FV->getSize ();
930  for (n = 0; n < len; n++) {
931  // check iteration voltages
932  nr_double_t v_abs = abs (VS->get (n) - VP->get (n));
933  nr_double_t v_rel = abs (VS->get (n));
934  if (v_abs >= vabstol + reltol * v_rel) return 0;
935  // check balanced currents
936  nr_complex_t il = IL->get (n);
937  nr_complex_t in = IN->get (n);
938  if (il != in) {
939  nr_double_t i_abs = abs (il + in);
940  nr_double_t i_rel = abs ((il + in) / (il - in));
941  if (i_abs >= iabstol && 2.0 * i_rel >= reltol) return 0;
942  }
943  }
944  return 1;
945 }
946 
947 // some definitions for the non-linear matrix filler
948 #undef G_
949 #undef C_
950 #define G_(r,c) (*jg) ((r)*nlfreqs+f,(c)*nlfreqs+f)
951 #define C_(r,c) (*jq) ((r)*nlfreqs+f,(c)*nlfreqs+f)
952 #undef FI_
953 #undef FQ_
954 #define FI_(r) (*ig) ((r)*nlfreqs+f)
955 #define FQ_(r) (*fq) ((r)*nlfreqs+f)
956 #define IR_(r) (*ir) ((r)*nlfreqs+f)
957 #define QR_(r) (*qr) ((r)*nlfreqs+f)
958 
959 /* This function fills in the matrix and vector entries for the
960  non-linear HB equations for a given frequency index. */
967  int f) {
968  // through each linear circuit
969  for (ptrlistiterator<circuit> it (nolcircuits); *it; ++it) {
970  circuit * cir = it.current ();
971  int s = cir->getSize ();
972  int nr, nc, r, c;
973 
974  for (r = 0; r < s; r++) {
975  if ((nr = cir->getNode(r)->getNode () - 1) < 0) continue;
976  // apply G- and C-matrix entries
977  for (c = 0; c < s; c++) {
978  if ((nc = cir->getNode(c)->getNode () - 1) < 0) continue;
979  G_(nr, nc) += cir->getY (r, c);
980  C_(nr, nc) += cir->getQV (r, c);
981  }
982  // apply I- and Q-vector entries
983  FI_(nr) -= cir->getI (r);
984  FQ_(nr) -= cir->getQ (r);
985  // ThinkME: positive or negative?
986  IR_(nr) += cir->getGV (r) + cir->getI (r);
987  QR_(nr) += cir->getCV (r) + cir->getQ (r);
988  }
989  }
990 }
991 
992 /* The function initializes the non-linear part of the HB. */
994  int N = nbanodes;
995 
996  // allocate matrices and vectors
997  if (FQ == NULL) {
998  FQ = new tvector<nr_complex_t> (N * nlfreqs);
999  }
1000  if (IG == NULL) {
1001  IG = new tvector<nr_complex_t> (N * nlfreqs);
1002  }
1003  if (IR == NULL) {
1004  IR = new tvector<nr_complex_t> (N * nlfreqs);
1005  }
1006  if (QR == NULL) {
1007  QR = new tvector<nr_complex_t> (N * nlfreqs);
1008  }
1009  if (JG == NULL) {
1010  JG = new tmatrix<nr_complex_t> (N * nlfreqs);
1011  }
1012  if (JQ == NULL) {
1013  JQ = new tmatrix<nr_complex_t> (N * nlfreqs);
1014  }
1015  if (JF == NULL) {
1016  JF = new tmatrix<nr_complex_t> (N * nlfreqs);
1017  }
1018 
1019  // voltage vector in frequency and time domain
1020  if (VS == NULL) {
1021  VS = new tvector<nr_complex_t> (N * nlfreqs);
1022  }
1023  if (vs == NULL) {
1024  vs = new tvector<nr_complex_t> (N * nlfreqs);
1025  }
1026  if (VP == NULL) {
1027  VP = new tvector<nr_complex_t> (N * nlfreqs);
1028  }
1029 
1030  // error vector
1031  if (FV == NULL) {
1032  FV = new tvector<nr_complex_t> (N * nlfreqs);
1033  }
1034  // right hand side vector
1035  if (RH == NULL) {
1036  RH = new tvector<nr_complex_t> (N * nlfreqs);
1037  }
1038 
1039  // linear and non-linear current vector
1040  if (IL == NULL) {
1041  IL = new tvector<nr_complex_t> (N * nlfreqs);
1042  }
1043  if (IN == NULL) {
1044  IN = new tvector<nr_complex_t> (N * nlfreqs);
1045  }
1046 
1047  // assign nodes
1049 
1050  // initialize circuits
1051  for (ptrlistiterator<circuit> it (nolcircuits); *it; ++it) {
1052  (*it)->initHB (nlfreqs);
1053  }
1054 }
1055 
1056 /* Saves the node voltages of the given circuit and for the given
1057  frequency entry into the circuit voltage vector. */
1059  int r, nr, s = cir->getSize ();
1060  for (r = 0; r < s; r++) {
1061  if ((nr = cir->getNode(r)->getNode () - 1) < 0) continue;
1062  // apply V-vector entries
1063  cir->setV (r, real (vs->get (nr * nlfreqs + f)));
1064  }
1065 }
1066 
1067 /* The function saves voltages into non-linear circuits, runs each
1068  non-linear components' HB calculator for each frequency and applies
1069  the matrix and vector entries appropriately. */
1071  // clear matrices and vectors before
1072  IG->set (0.0);
1073  FQ->set (0.0);
1074  IR->set (0.0);
1075  QR->set (0.0);
1076  JG->set (0.0);
1077  JQ->set (0.0);
1078  // through each frequency
1079  for (int f = 0; f < nlfreqs; f++) {
1080  // calculate components' HB matrices and vector for the given frequency
1081  for (ptrlistiterator<circuit> it (nolcircuits); *it; ++it) {
1082  saveNodeVoltages (*it, f); // node voltages
1083  (*it)->calcHB (f); // HB calculator
1084  }
1085  // fill in all matrix entries for the given frequency
1086  fillMatrixNonLinear (JG, JQ, IG, FQ, IR, QR, f);
1087  }
1088 }
1089 
1090 /* The following function transforms a vector using a Fast Fourier
1091  Transformation from the time domain to the frequency domain. */
1093  int i, k, r;
1094  int n = nlfreqs;
1095  int nd = dfreqs.getSize ();
1096  int nodes = V->getSize () / n;
1097  nr_double_t * d = (nr_double_t *) V->getData ();
1098 
1099  if (nd == 1) {
1100  // for each node a single 1d-FFT
1101  for (k = i = 0; i < nodes; i++, k += 2 * n) {
1102  nr_double_t * dst = &d[k];
1103  _fft_1d (dst, n, isign);
1104  if (isign > 0) for (r = 0; r < 2 * n; r++) *dst++ /= n;
1105  }
1106  }
1107  else {
1108  // for each node a single nd-FFT
1109  for (k = i = 0; i < nodes; i++, k += 2 * n) {
1110  nr_double_t * dst = &d[k];
1111  _fft_nd (dst, ndfreqs, nd, isign);
1112  if (isign > 0) for (r = 0; r < 2 * n; r++) *dst++ /= ndfreqs[0];
1113  }
1114  }
1115 }
1116 
1117 /* The following function transforms a vector using an Inverse Fast
1118  Fourier Transformation from the frequency domain to the domain
1119  time. */
1121  VectorFFT (V, -isign);
1122 }
1123 
1124 /* The following function transforms a matrix using a Fast Fourier
1125  Transformation from the time domain to the frequency domain. */
1127 
1128 #if THE_SLO_ALGO
1129  // each column FFT
1130  for (int c = 0; c < M->getCols (); c++) {
1132  VectorFFT (&V);
1133  M->setCol (c, V);
1134  }
1135 
1136  // each row IFFT
1137  for (int r = 0; r < M->getRows (); r++) {
1139  VectorIFFT (&V);
1140  M->setRow (r, V);
1141  }
1142 #else
1143  int c, r, nc, nr;
1144  // for each non-linear node block
1145  for (nc = c = 0; c < nbanodes; c++, nc += nlfreqs) {
1146  for (nr = r = 0; r < nbanodes; r++, nr += nlfreqs) {
1148  int fr, fc, fi;
1149  // transform the sub-diagonal only
1150  for (fc = 0; fc < nlfreqs; fc++) V (fc) = M->get (nr + fc, nc + fc);
1151  VectorFFT (&V);
1152  // fill in resulting sub-matrix for the node
1153  for (fc = 0; fc < nlfreqs; fc++) {
1154  for (fi = nlfreqs - 1 - fc, fr = 0; fr < nlfreqs; fr++) {
1155  if (++fi >= nlfreqs) fi = 0;
1156  M->set (nr + fr, nc + fc, V (fi));
1157  }
1158  }
1159  }
1160  }
1161 #endif
1162 }
1163 
1164 /* This function solves the actual HB equation in the frequency domain.
1165  F(V) = IC + [YV] * VS + j[O] * FQ + IG -> 0
1166  Care must be taken with indexing here: In the frequency domain only
1167  real positive frequencies are used and computed, but in the time
1168  domain we have more values because of the periodic continuation in
1169  the frequency domain.
1170  RHS = j[O] * CV + GV - (IC + j[O] * FQ + IG)
1171  Also the right hand side of the equation system for the new voltage
1172  vector is computed here. */
1173 void hbsolver::solveHB (void) {
1174  // for each non-linear node
1175  for (int r = 0; r < nbanodes * nlfreqs; ) {
1176  // for each frequency
1177  for (int f = 0; f < nlfreqs; f++, r++) {
1178  nr_complex_t il = 0.0, in = 0.0, ir = 0.0;
1179  // constant current vector due to sources
1180  il += IC->get (r);
1181  // part 1 of right hand side vector
1182  ir -= il;
1183  // transadmittance matrix multiplied by voltage vector
1184  for (int c = 0; c < nbanodes * nlfreqs; c++) {
1185  il += YV_(r, c) * VS_(c);
1186  }
1187  // charge vector
1188  in += OM_(f) * FQ->get (r);
1189  // current vector
1190  in += IG->get (r);
1191  // part 2, 3 and 4 of right hand side vector
1192  ir += IR->get (r);
1193  ir += OM_(f) * QR->get (r);
1194  // put values into result vectors
1195  RH->set (r, ir);
1196  FV->set (r, il + in);
1197  IL->set (r, il);
1198  IN->set (r, in);
1199  }
1200  }
1201 }
1202 
1203 /* The function calculates the full Jacobian JF = [YV] + j[O] * JQ + JG */
1205  int c, r, fc, fr, rt, ct;
1206  /* add admittances of capacitance matrix JQ and non-linear
1207  admittances matrix JG into complete Jacobian JF */
1208  for (c = 0; c < nbanodes; c++) {
1209  ct = c * nlfreqs;
1210  for (fc = 0; fc < nlfreqs; fc++, ct++) {
1211  for (r = 0; r < nbanodes; r++) {
1212  rt = r * nlfreqs;
1213  for (fr = 0; fr < nlfreqs; fr++, rt++) {
1214  JF_(rt, ct) = JG->get (rt, ct) + JQ->get (rt, ct) * OM_(fr);
1215  }
1216  }
1217  }
1218  }
1219  *JF += *YV; // add linear admittance matrix
1220 }
1221 
1222 /* The function expands the given vector in the frequency domain to
1223  make it a real valued signal in the time domain. */
1225  int nodes) {
1226  tvector<nr_complex_t> res (nodes * nlfreqs);
1227  int r, ff, rf, rt;
1228  for (r = 0; r < nodes; r++) {
1229  rt = r * nlfreqs;
1230  rf = r * lnfreqs;
1231  // copy first part of vector
1232  for (ff = 0; ff < lnfreqs; ff++, rf++, rt++) {
1233  res (rt) = V (rf);
1234  }
1235  // continue vector conjugated
1236  for (rf -= 2; ff < nlfreqs; ff++, rf--, rt++) {
1237  res (rt) = conj (V (rf));
1238  }
1239  }
1240  return res;
1241 }
1242 
1243 /* The function expands the given matrix in the frequency domain to
1244  make it a real valued signal in the time domain. */
1246  int nodes) {
1247  tmatrix<nr_complex_t> res (nodes * nlfreqs);
1248  int r, c, rf, rt, cf, ct, ff;
1249  for (r = 0; r < nodes; r++) {
1250  for (c = 0; c < nodes; c++) {
1251  rf = r * lnfreqs;
1252  rt = r * nlfreqs;
1253  cf = c * lnfreqs;
1254  ct = c * nlfreqs;
1255  // copy first part of diagonal
1256  for (ff = 0; ff < lnfreqs; ff++, cf++, ct++, rf++, rt++) {
1257  res (rt, ct) = M (rf, cf);
1258  }
1259  // continue diagonal conjugated
1260  for (cf -= 2, rf -= 2; ff < nlfreqs; ff++, cf--, ct++, rf--, rt++) {
1261  res (rt, ct) = conj (M (rf, cf));
1262  }
1263  }
1264  }
1265  return res;
1266 }
1267 
1268 /* This function solves the equation system
1269  JF * VS(n+1) = JF * VS(n) - FV
1270  in order to obtains a new voltage vector in the frequency domain. */
1272  // save previous iteration voltage
1273  *VP = *VS;
1274 
1275  // setup equation system
1277  try_running () {
1278  // use LU decomposition for solving
1280  eqns.passEquationSys (JF, VS, RH);
1281  eqns.solve ();
1282  }
1283  // appropriate exception handling
1284  catch_exception () {
1285  case EXCEPTION_PIVOT:
1286  default:
1287  logprint (LOG_ERROR, "WARNING: %s: during NR iteration\n", getName ());
1288  estack.print ();
1289  }
1290 
1291  // save new voltages in time domain vector
1292  *vs = *VS;
1293 }
1294 
1295 /* The following function extends the existing linear MNA matrix to
1296  contain the additional rows and columns for the excitation voltage
1297  sources. */
1299  int nodes) {
1300  int no = M.getCols ();
1301  tmatrix<nr_complex_t> res (no + nodes * lnfreqs);
1302  // copy the existing part
1303  for (int r = 0; r < no; r++) {
1304  for (int c = 0; c < no; c++) {
1305  res (r, c) = M (r, c);
1306  }
1307  }
1308  return res;
1309 }
1310 
1311 /* The function fills in the missing MNA entries for the excitation
1312  voltage sources into the extended rows and columns as well as the
1313  actual voltage values into the right hand side vector. */
1315  tvector<nr_complex_t> * I) {
1316  // through each excitation source
1317  int of = lnfreqs * (nlnvsrcs + nnanodes);
1318  int sc = of;
1319 
1320  for (ptrlistiterator<circuit> it (excitations); *it; ++it) {
1321  circuit * vs = it.current ();
1322  // get positive and negative node
1323  int pnode = vs->getNode(NODE_1)->getNode ();
1324  int nnode = vs->getNode(NODE_2)->getNode ();
1325  for (int f = 0; f < lnfreqs; f++, sc++) { // for each frequency
1326  // fill right hand side vector
1327  nr_double_t freq = rfreqs (f);
1328  vs->calcHB (freq);
1329  I_(sc) = vs->getE (VSRC_1);
1330  // fill MNA entries
1331  int pn = (pnode - 1) * lnfreqs + f;
1332  int nn = (nnode - 1) * lnfreqs + f;
1333  if (pnode) {
1334  Y_(pn, sc) = +1.0;
1335  Y_(sc, pn) = +1.0;
1336  }
1337  if (nnode) {
1338  Y_(nn, sc) = -1.0;
1339  Y_(sc, nn) = -1.0;
1340  }
1341  }
1342  }
1343 }
1344 
1345 /* The function calculates and saves the final solution. */
1347 
1348  // extend the linear MNA matrix
1350 
1351  int S = NA->getCols ();
1352  int N = nnanodes * lnfreqs;
1353 
1354  // right hand side vector
1356  // temporary solution
1358  // final solution
1359  x = new tvector<nr_complex_t> (N);
1360 
1361  // fill in missing MNA entries
1363 
1364  // put currents through balanced nodes into right hand side
1365  for (int n = 0; n < nbanodes; n++) {
1366  for (int f = 0; f < lnfreqs; f++) {
1367  nr_complex_t i = IL->get (n * nlfreqs + f);
1368  if (f != 0 && f != lnfreqs - 1) i *= 2;
1369  I_(n * lnfreqs + f) = i;
1370  }
1371  }
1372 
1373  // use QR decomposition for the final solution
1374  try_running () {
1377  eqns.passEquationSys (NA, V, I);
1378  eqns.solve ();
1379  }
1380  // appropriate exception handling
1381  catch_exception () {
1382  case EXCEPTION_PIVOT:
1383  default:
1384  logprint (LOG_ERROR, "WARNING: %s: during final AC analysis\n",
1385  getName ());
1386  estack.print ();
1387  }
1388  for (int i = 0; i < N; i++) x->set (i, V_(i));
1389 }
1390 
1391 // Saves simulation results.
1393  vector * f;
1394  // add current frequency to the dependency of the output dataset
1395  if ((f = data->findDependency ("hbfrequency")) == NULL) {
1396  f = new vector ("hbfrequency");
1397  data->addDependency (f);
1398  }
1399  // save frequency vector
1400  if (runs == 1) {
1401  for (int i = 0; i < lnfreqs; i++) f->add (rfreqs (i));
1402  }
1403  // save node voltage vectors
1404  int nanode = 0;
1405  for (strlistiterator it (nanodes); *it; ++it, nanode++) {
1406  int l = strlen (*it);
1407  char * n = (char *) malloc (l + 4);
1408  sprintf (n, "%s.Vb", *it);
1409  for (int i = 0; i < lnfreqs; i++) {
1410  saveVariable (n, x->get (i + nanode * lnfreqs), f);
1411  }
1412  }
1413 }
1414 
1415 // properties
1416 PROP_REQ [] = {
1417  { "n", PROP_INT, { 1, PROP_NO_STR }, PROP_MIN_VAL (1) },
1418  PROP_NO_PROP };
1419 PROP_OPT [] = {
1420  { "f", PROP_REAL, { 1e9, PROP_NO_STR }, PROP_POS_RANGEX },
1421  { "iabstol", PROP_REAL, { 1e-12, PROP_NO_STR }, PROP_RNG_X01I },
1422  { "vabstol", PROP_REAL, { 1e-6, PROP_NO_STR }, PROP_RNG_X01I },
1423  { "reltol", PROP_REAL, { 1e-3, PROP_NO_STR }, PROP_RNG_X01I },
1424  { "MaxIter", PROP_INT, { 150, PROP_NO_STR }, PROP_RNGII (2, 10000) },
1425  PROP_NO_PROP };
1426 struct define_t hbsolver::anadef =
1428 
1429 } // namespace qucs
1430 
1431 /* TODO list for HB Solver:
1432  - Take care about nodes with non-linear components only.
1433  - AC Power Sources (extra Z and open loop voltage).
1434  - Current sources.
1435  - Balancing of multi-dimensional non-linear networks.
1436  - Sources directly connected to non-linear components and no other
1437  linear component (insert short).
1438  - Bug: With capacitors at hand there is voltage convergence but no
1439  current balancing.
1440  - Output currents through voltage sources.
1441  */
void add(char *)
Definition: strlist.cpp:65
type_t * current(void)
Definition: ptrlist.cpp:210
std::complex< nr_double_t > nr_complex_t
Definition: complex.h:31
void _fft_1d(nr_double_t *, int, int isign=1)
Definition: fourier.cpp:50
exceptionstack estack
tmatrix< nr_complex_t > * YV
Definition: hbsolver.h:103
l
Definition: parse_vcd.y:213
circuit * getRoot(void)
Definition: net.h:46
net * subnet
Definition: analysis.h:273
bool isExcitation(circuit *)
Definition: hbsolver.cpp:324
#define NODE_2
Definition: circuit.h:35
matrix real(matrix a)
Real part matrix.
Definition: matrix.cpp:568
matrix abs(matrix a)
Computes magnitude of each matrix element.
Definition: matrix.cpp:531
#define PROP_RNGII(f, t)
Definition: netdefs.h:138
#define N(n)
Definition: equation.cpp:58
#define I_(r)
Definition: hbsolver.cpp:685
tmatrix< nr_complex_t > * Y
Definition: hbsolver.h:99
strlist * banodes
Definition: hbsolver.h:94
void set(int, nr_type_t)
Definition: tvector.cpp:108
void saveVariable(const char *, nr_complex_t, qucs::vector *)
Save variable into analysis dataset.
Definition: analysis.cpp:151
int checkBalance(void)
Definition: hbsolver.cpp:925
#define PROP_DEF
Definition: netdefs.h:189
int assignVoltageSources(ptrlist< circuit >)
Definition: hbsolver.cpp:513
PROP_OPT[]
Definition: acsolver.cpp:232
void setRow(int, tvector< nr_type_t >)
Definition: tmatrix.cpp:144
nr_double_t getPropertyDouble(const char *)
Definition: object.cpp:176
void initDC(void)
Definition: hbsolver.cpp:316
#define A_(r, c)
Definition: hbsolver.cpp:687
#define C_(r, c)
Definition: hbsolver.cpp:951
void set(int, int, nr_type_t)
Definition: tmatrix.cpp:120
void VectorFFT(tvector< nr_complex_t > *, int isign=1)
Definition: hbsolver.cpp:1092
#define PROP_REAL
Definition: netdefs.h:174
int calcOrder(int)
Definition: hbsolver.cpp:365
tvector< nr_complex_t > * FQ
Definition: hbsolver.h:110
void MatrixFFT(tmatrix< nr_complex_t > *)
Definition: hbsolver.cpp:1126
int contains(char *)
Definition: strlist.cpp:107
tvector< nr_complex_t > * RH
Definition: hbsolver.h:119
void loadMatrices(void)
Definition: hbsolver.cpp:1070
#define PROP_NO_PROP
Definition: netdefs.h:122
strlist * exnodes
Definition: hbsolver.h:94
void collectFrequencies(void)
Definition: hbsolver.cpp:374
tmatrix< nr_complex_t > * NA
Definition: hbsolver.h:104
tvector< nr_complex_t > * IR
Definition: hbsolver.h:117
tvector< nr_complex_t > * IL
Definition: hbsolver.h:114
int getPropertyInteger(const char *)
Definition: object.cpp:198
void calcJacobian(void)
Definition: hbsolver.cpp:1204
void calcConstantCurrent(void)
Definition: hbsolver.cpp:871
tvector< nr_complex_t > * vs
Definition: hbsolver.h:125
#define PROP_NO_STR
Definition: netdefs.h:125
int getVoltageSource(void)
Definition: circuit.h:187
void passEquationSys(tmatrix< nr_type_t > *, tvector< nr_type_t > *, tvector< nr_type_t > *)
Definition: eqnsys.cpp:99
tvector< nr_complex_t > * FV
Definition: hbsolver.h:113
#define FI_(r)
Definition: hbsolver.cpp:954
#define catch_exception()
void setVoltageSource(int s)
Definition: circuit.h:186
n
Definition: parse_citi.y:147
#define S(con)
Definition: property.cpp:158
#define PROP_LINEAR
Definition: netdefs.h:120
tvector< nr_complex_t > * IC
Definition: hbsolver.h:122
r
Definition: parse_mdl.y:515
#define PROP_INT
Definition: netdefs.h:173
#define Y_(r, c)
Definition: hbsolver.cpp:690
#define ZCU_(r, c)
Definition: hbsolver.cpp:694
void saveResults(void)
Definition: hbsolver.cpp:1392
#define PROP_RNG_X01I
Definition: netdefs.h:137
int getSize(void)
Get the number of ports the circuit element has.
Definition: circuit.h:143
#define throw_exception(e)
#define PROP_ACTION
Definition: netdefs.h:115
#define try_running()
void VectorIFFT(tvector< nr_complex_t > *, int isign=1)
Definition: hbsolver.cpp:1120
tvector< nr_type_t > getRow(int)
Definition: tmatrix.cpp:133
void invertMatrix(tmatrix< nr_complex_t > *, tmatrix< nr_complex_t > *)
Definition: hbsolver.cpp:649
#define VSRC_1
Definition: circuit.h:40
i
Definition: parse_mdl.y:516
hbsolver(char *)
Definition: hbsolver.cpp:70
int getVoltageSources(void)
Definition: circuit.cpp:602
void expandFrequencies(nr_double_t, int)
Definition: hbsolver.cpp:332
void getNodeLists(void)
Definition: hbsolver.cpp:471
tmatrix< nr_complex_t > * JF
Definition: hbsolver.h:108
ptrlist< circuit > nolcircuits
Definition: hbsolver.h:96
void add(nr_type_t)
Definition: tvector.cpp:127
nr_complex_t getC(int, int)
Definition: circuit.cpp:355
int getNode(void)
Definition: node.cpp:71
virtual void calcHB(nr_double_t)
Definition: circuit.h:125
int assignNodes(ptrlist< circuit >, strlist *, int offset=0)
Definition: hbsolver.cpp:526
#define ZVL_(r, c)
Definition: hbsolver.cpp:693
tmatrix< nr_complex_t > * JG
Definition: hbsolver.h:107
strlist * circuitNodes(ptrlist< circuit >)
Definition: hbsolver.cpp:456
base class for qucs circuit elements.
Definition: circuit.h:92
int getType(void)
Definition: circuit.h:137
static void calc(hbsolver *)
Definition: hbsolver.cpp:298
strlist * lnnodes
Definition: hbsolver.h:94
#define IR_(r)
Definition: hbsolver.cpp:956
nr_complex_t getI(int)
Definition: circuit.cpp:391
tvector< nr_double_t > posfreqs
Definition: hbsolver.h:89
void setCol(int, tvector< nr_type_t >)
Definition: tmatrix.cpp:164
tvector< nr_double_t > dfreqs
Definition: hbsolver.h:92
nr_type_t get(int)
Definition: tvector.cpp:101
void add(nr_complex_t)
Definition: vector.cpp:151
net * getNet(void)
Definition: circuit.h:173
#define V(con)
Definition: evaluate.cpp:63
void solveVoltages(void)
Definition: hbsolver.cpp:1271
void prepareLinear(void)
Definition: hbsolver.cpp:546
void fillMatrixLinearExtended(tmatrix< nr_complex_t > *, tvector< nr_complex_t > *)
Definition: hbsolver.cpp:1314
void splitCircuits(void)
Definition: hbsolver.cpp:437
int solve(void)
placehoder for solution function
Definition: hbsolver.cpp:145
The analysis class header file.
#define VS_(r)
Definition: hbsolver.cpp:140
x
Definition: parse_mdl.y:498
void createMatrixLinearY(void)
Definition: hbsolver.cpp:711
int length(void)
Definition: strlist.cpp:100
void prepareNonLinear(void)
Definition: hbsolver.cpp:993
nr_complex_t getE(int)
Definition: circuit.cpp:379
dataset * data
Definition: analysis.h:274
void print(const char *prefix=NULL)
#define M_PI
Archimedes' constant ( )
Definition: consts.h:47
void solve(void)
Definition: eqnsys.cpp:124
void setAlgo(int a)
Definition: eqnsys.h:69
#define FQ_(r)
Definition: hbsolver.cpp:955
tvector< nr_double_t > rfreqs
Definition: hbsolver.h:90
void saveNodeVoltages(circuit *, int)
Definition: hbsolver.cpp:1058
tmatrix< nr_complex_t > * Z
Definition: hbsolver.h:101
The circuit class header file.
#define YV_(r, c)
Definition: hbsolver.cpp:697
nr_complex_t getD(int, int)
Definition: circuit.cpp:367
tvector< nr_complex_t > * VS
Definition: hbsolver.h:111
tvector< nr_complex_t > * x
Definition: hbsolver.h:124
nodes
nr_complex_t getGV(int)
Definition: circuit.cpp:493
type
Definition: parse_vcd.y:164
void _fft_nd(nr_double_t *, int[], int, int isign=1)
Definition: fourier.cpp:244
#define PROP_MIN_VAL(k)
Definition: netdefs.h:133
tvector< nr_complex_t > * VP
Definition: hbsolver.h:112
void initHB(void)
Definition: hbsolver.cpp:307
#define JF_(r, c)
Definition: hbsolver.cpp:699
strlist * nanodes
Definition: hbsolver.h:94
nr_double_t frequency
Definition: hbsolver.h:93
void print(bool realonly=false)
v
Definition: parse_zvr.y:141
tmatrix< nr_complex_t > expandMatrix(tmatrix< nr_complex_t >, int)
Definition: hbsolver.cpp:1245
#define ZVU_(r, c)
Definition: hbsolver.cpp:692
char * toString(const char *concat=" ")
Definition: strlist.cpp:177
tvector< nr_complex_t > * QR
Definition: hbsolver.h:118
eqns
nr_complex_t getB(int, int)
Definition: circuit.cpp:343
tvector< nr_complex_t > * IS
Definition: hbsolver.h:123
nr_type_t get(int, int)
Definition: tmatrix.cpp:113
void fillMatrixLinearA(tmatrix< nr_complex_t > *, int)
Definition: hbsolver.cpp:597
int * ndfreqs
Definition: hbsolver.h:91
#define PROP_POS_RANGEX
Definition: netdefs.h:131
tvector< nr_type_t > getCol(int)
Definition: tmatrix.cpp:153
#define A(a)
Definition: eqndefined.cpp:72
virtual void initHB(void)
Definition: circuit.h:124
#define OM_(r)
Definition: hbsolver.cpp:141
nr_complex_t getQV(int, int)
Definition: circuit.cpp:481
void fillMatrixNonLinear(tmatrix< nr_complex_t > *, tmatrix< nr_complex_t > *, tvector< nr_complex_t > *, tvector< nr_complex_t > *, tvector< nr_complex_t > *, tvector< nr_complex_t > *, int)
Definition: hbsolver.cpp:961
#define V_(r)
Definition: hbsolver.cpp:684
tvector< nr_complex_t > expandVector(tvector< nr_complex_t >, int)
Definition: hbsolver.cpp:1224
#define B_(r, c)
Definition: hbsolver.cpp:591
#define NODE_1
Definition: circuit.h:34
void print(void)
void setText(const char *,...)
Definition: exception.cpp:67
tvector< nr_complex_t > * IG
Definition: hbsolver.h:109
tvector< nr_complex_t > * IN
Definition: hbsolver.h:115
#define LOG_ERROR
Definition: logging.h:28
void finalSolution(void)
Definition: hbsolver.cpp:1346
char * getName(void)
Definition: object.cpp:84
char * get(int)
Definition: strlist.cpp:129
nr_complex_t excitationZ(tvector< nr_complex_t > *, circuit *, int)
Definition: hbsolver.cpp:857
node * getNode(int)
Definition: circuit.cpp:307
#define G_(r, c)
Definition: hbsolver.cpp:950
matrix conj(matrix a)
Conjugate complex matrix.
Definition: matrix.cpp:505
ptrlist< circuit > excitations
Definition: hbsolver.h:95
#define LOG_STATUS
Definition: logging.h:29
tmatrix< nr_complex_t > * JQ
Definition: hbsolver.h:106
tvector< nr_double_t > negfreqs
Definition: hbsolver.h:88
PROP_REQ[]
Definition: acsolver.cpp:229
void clear(void)
Definition: tvector.cpp:164
#define ZCL_(r, c)
Definition: hbsolver.cpp:695
nr_complex_t getQ(int)
Definition: circuit.cpp:413
nr_complex_t getY(int, int)
Definition: circuit.cpp:446
void createMatrixLinearA(void)
Definition: hbsolver.cpp:563
#define M(con)
Definition: evaluate.cpp:64
void logprint(int level, const char *format,...)
Definition: logging.c:37
strlist * nlnodes
Definition: hbsolver.h:94
tmatrix< nr_complex_t > extendMatrixLinear(tmatrix< nr_complex_t >, int)
Definition: hbsolver.cpp:1298
int getCols(void)
Definition: tmatrix.h:60
ptrlist< circuit > lincircuits
Definition: hbsolver.h:97
void setNode(int)
Definition: node.cpp:66
tmatrix< nr_complex_t > * A
Definition: hbsolver.h:100
int getSize(void)
Definition: tvector.h:82
#define D_(r, c)
Definition: hbsolver.cpp:593
void solveHB(void)
Definition: hbsolver.cpp:1173
int contains(nr_type_t, nr_double_t eps=NR_EPSI)
Definition: tvector.cpp:171
nr_type_t * getData(void)
Definition: tvector.h:83
tvector< nr_complex_t > * OM
Definition: hbsolver.h:120
void setV(int, nr_complex_t)
Definition: circuit.cpp:440
nr_complex_t getCV(int)
Definition: circuit.cpp:505
int getRows(void)
Definition: tmatrix.h:61
#define QR_(r)
Definition: hbsolver.cpp:957