Qucs-core  0.0.18
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
eqnsys.cpp
Go to the documentation of this file.
1 /*
2  * eqnsys.cpp - equation system solver class implementation
3  *
4  * Copyright (C) 2004, 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 // the types required for qucs library files are defined
26 // in qucs_typedefs.h, created by configure from
27 // qucs_typedefs.h.in
28 #include "qucs_typedefs.h"
29 
30 #include <assert.h>
31 #include <time.h>
32 #include <cmath>
33 #include <float.h>
34 
35 #include "compat.h"
36 #include "logging.h"
37 #include "precision.h"
38 #include "complex.h"
39 #include "tmatrix.h"
40 #include "eqnsys.h"
41 #include "exception.h"
42 #include "exceptionstack.h"
43 
45 #define Swap(type,a,b) { type t; t = a; a = b; b = t; }
46 
47 namespace qucs {
48 
50 template <class nr_type_t>
52  A = V = NULL;
53  B = X = NULL;
54  S = E = NULL;
55  T = R = NULL;
56  nPvt = NULL;
57  cMap = rMap = NULL;
58  update = 1;
59  pivoting = PIVOT_PARTIAL;
60  N = 0;
61 }
62 
64 template <class nr_type_t>
66  if (R != NULL) delete R;
67  if (T != NULL) delete T;
68  if (B != NULL) delete B;
69  if (S != NULL) delete S;
70  if (E != NULL) delete E;
71  if (V != NULL) delete V;
72  if (rMap != NULL) delete[] rMap;
73  if (cMap != NULL) delete[] cMap;
74  if (nPvt != NULL) delete[] nPvt;
75 }
76 
79 template <class nr_type_t>
81  A = e.A;
82  V = NULL;
83  S = E = NULL;
84  T = R = NULL;
85  B = e.B ? new tvector<nr_type_t> (*(e.B)) : NULL;
86  cMap = rMap = NULL;
87  nPvt = NULL;
88  update = 1;
89  X = e.X;
90  N = 0;
91 }
92 
98 template <class nr_type_t>
100  tvector<nr_type_t> * refX,
101  tvector<nr_type_t> * nB) {
102  if (nA != NULL) {
103  A = nA;
104  update = 1;
105  if (N != A->getCols ()) {
106  N = A->getCols ();
107  if (cMap) delete[] cMap; cMap = new int[N];
108  if (rMap) delete[] rMap; rMap = new int[N];
109  if (nPvt) delete[] nPvt; nPvt = new nr_double_t[N];
110  }
111  }
112  else {
113  update = 0;
114  }
115  if (B != NULL) delete B;
116  B = new tvector<nr_type_t> (*nB);
117  X = refX;
118 }
119 
123 template <class nr_type_t>
125 #if DEBUG && 0
126  time_t t = time (NULL);
127 #endif
128  switch (algo) {
129  case ALGO_INVERSE:
130  solve_inverse ();
131  break;
132  case ALGO_GAUSS:
133  solve_gauss ();
134  break;
135  case ALGO_GAUSS_JORDAN:
136  solve_gauss_jordan ();
137  break;
139  solve_lu_crout ();
140  break;
142  solve_lu_doolittle ();
143  break;
145  factorize_lu_crout ();
146  break;
148  factorize_lu_doolittle ();
149  break;
151  substitute_lu_crout ();
152  break;
154  substitute_lu_doolittle ();
155  break;
156  case ALGO_JACOBI: case ALGO_GAUSS_SEIDEL:
157  solve_iterative ();
158  break;
159  case ALGO_SOR:
160  solve_sor ();
161  break;
163  solve_qr ();
164  break;
166  solve_qr_ls ();
167  break;
169  solve_svd ();
170  break;
172  solve_qrh ();
173  break;
174  }
175 #if DEBUG && 0
176  logprint (LOG_STATUS, "NOTIFY: %dx%d eqnsys solved in %ld seconds\n",
177  A->getRows (), A->getCols (), time (NULL) - t);
178 #endif
179 }
180 
182 template <class nr_type_t>
184  *X = inverse (*A) * *B;
185 }
186 
187 #define A_ (*A)
188 #define X_ (*X)
189 #define B_ (*B)
190 
193 template <class nr_type_t>
195  nr_double_t MaxPivot;
196  nr_type_t f;
197  int i, c, r, pivot;
198 
199  // triangulate the matrix
200  for (i = 0; i < N; i++) {
201  // find maximum column value for pivoting
202  for (MaxPivot = 0, pivot = r = i; r < N; r++) {
203  if (abs (A_(r, i)) > MaxPivot) {
204  MaxPivot = abs (A_(r, i));
205  pivot = r;
206  }
207  }
208  // exchange rows if necessary
209  assert (MaxPivot != 0);
210  if (i != pivot) {
211  A->exchangeRows (i, pivot);
212  B->exchangeRows (i, pivot);
213  }
214  // compute new rows and columns
215  for (r = i + 1; r < N; r++) {
216  f = A_(r, i) / A_(i, i);
217  for (c = i + 1; c < N; c++) A_(r, c) -= f * A_(i, c);
218  B_(r) -= f * B_(i);
219  }
220  }
221 
222  // backward substitution
223  for (i = N - 1; i >= 0; i--) {
224  f = B_(i);
225  for (c = i + 1; c < N; c++) f -= A_(i, c) * X_(c);
226  X_(i) = f / A_(i, i);
227  }
228 }
229 
233 template <class nr_type_t>
235  nr_double_t MaxPivot;
236  nr_type_t f;
237  int i, c, r, pivot;
238 
239  // create the eye matrix
240  for (i = 0; i < N; i++) {
241  // find maximum column value for pivoting
242  for (MaxPivot = 0, pivot = r = i; r < N; r++) {
243  if (abs (A_(r, i)) > MaxPivot) {
244  MaxPivot = abs (A_(r, i));
245  pivot = r;
246  }
247  }
248  // exchange rows if necessary
249  assert (MaxPivot != 0);
250  if (i != pivot) {
251  A->exchangeRows (i, pivot);
252  B->exchangeRows (i, pivot);
253  }
254 
255  // compute current row
256  f = A_(i, i);
257  for (c = i + 1; c < N; c++) A_(i, c) /= f;
258  B_(i) /= f;
259 
260  // compute new rows and columns
261  for (r = 0; r < N; r++) {
262  if (r != i) {
263  f = A_(r, i);
264  for (c = i + 1; c < N; c++) A_(r, c) -= f * A_(i, c);
265  B_(r) -= f * B_(i);
266  }
267  }
268  }
269 
270  // right hand side is now the solution
271  *X = *B;
272 }
273 
274 #define LU_FAILURE 0
275 #define VIRTUAL_RES(txt,i) { \
276  qucs::exception * e = new qucs::exception (EXCEPTION_SINGULAR); \
277  e->setText (txt); \
278  e->setData (rMap[i]); \
279  A_(i, i) = NR_TINY; /* virtual resistance to ground */ \
280  throw_exception (e); }
281 
287 template <class nr_type_t>
289 
290  // skip decomposition if requested
291  if (update) {
292  // perform LU composition
293  factorize_lu_crout ();
294  }
295 
296  // finally solve the equation system
297  substitute_lu_crout ();
298 }
299 
301 template <class nr_type_t>
303 
304  // skip decomposition if requested
305  if (update) {
306  // perform LU composition
307  factorize_lu_doolittle ();
308  }
309 
310  // finally solve the equation system
311  substitute_lu_doolittle ();
312 }
313 
318 template <class nr_type_t>
320  nr_double_t d, MaxPivot;
321  nr_type_t f;
322  int k, c, r, pivot;
323 
324  // initialize pivot exchange table
325  for (r = 0; r < N; r++) {
326  for (MaxPivot = 0, c = 0; c < N; c++)
327  if ((d = abs (A_(r, c))) > MaxPivot)
328  MaxPivot = d;
329  if (MaxPivot <= 0) MaxPivot = NR_TINY;
330  nPvt[r] = 1 / MaxPivot;
331  rMap[r] = r;
332  }
333 
334  // decompose the matrix into L (lower) and U (upper) matrix
335  for (c = 0; c < N; c++) {
336  // upper matrix entries
337  for (r = 0; r < c; r++) {
338  f = A_(r, c);
339  for (k = 0; k < r; k++) f -= A_(r, k) * A_(k, c);
340  A_(r, c) = f / A_(r, r);
341  }
342  // lower matrix entries
343  for (MaxPivot = 0, pivot = r; r < N; r++) {
344  f = A_(r, c);
345  for (k = 0; k < c; k++) f -= A_(r, k) * A_(k, c);
346  A_(r, c) = f;
347  // larger pivot ?
348  if ((d = nPvt[r] * abs (f)) > MaxPivot) {
349  MaxPivot = d;
350  pivot = r;
351  }
352  }
353 
354  // check pivot element and throw appropriate exception
355  if (MaxPivot <= 0) {
356 #if LU_FAILURE
358  e->setText ("no pivot != 0 found during Crout LU decomposition");
359  e->setData (c);
360  throw_exception (e);
361  goto fail;
362 #else /* insert virtual resistance */
363  VIRTUAL_RES ("no pivot != 0 found during Crout LU decomposition", c);
364 #endif
365  }
366 
367  // swap matrix rows if necessary and remember that step in the
368  // exchange table
369  if (c != pivot) {
370  A->exchangeRows (c, pivot);
371  Swap (int, rMap[c], rMap[pivot]);
372  Swap (nr_double_t, nPvt[c], nPvt[pivot]);
373  }
374  }
375 #if LU_FAILURE
376  fail:
377 #endif
378 }
379 
384 template <class nr_type_t>
386  nr_double_t d, MaxPivot;
387  nr_type_t f;
388  int k, c, r, pivot;
389 
390  // initialize pivot exchange table
391  for (r = 0; r < N; r++) {
392  for (MaxPivot = 0, c = 0; c < N; c++)
393  if ((d = abs (A_(r, c))) > MaxPivot)
394  MaxPivot = d;
395  if (MaxPivot <= 0) MaxPivot = NR_TINY;
396  nPvt[r] = 1 / MaxPivot;
397  rMap[r] = r;
398  }
399 
400  // decompose the matrix into L (lower) and U (upper) matrix
401  for (c = 0; c < N; c++) {
402  // upper matrix entries
403  for (r = 0; r < c; r++) {
404  f = A_(r, c);
405  for (k = 0; k < r; k++) f -= A_(r, k) * A_(k, c);
406  A_(r, c) = f;
407  }
408  // lower matrix entries
409  for (MaxPivot = 0, pivot = r; r < N; r++) {
410  f = A_(r, c);
411  for (k = 0; k < c; k++) f -= A_(r, k) * A_(k, c);
412  A_(r, c) = f;
413  // larger pivot ?
414  if ((d = nPvt[r] * abs (f)) > MaxPivot) {
415  MaxPivot = d;
416  pivot = r;
417  }
418  }
419 
420  // check pivot element and throw appropriate exception
421  if (MaxPivot <= 0) {
422 #if LU_FAILURE
424  e->setText ("no pivot != 0 found during Doolittle LU decomposition");
425  e->setData (c);
426  throw_exception (e);
427  goto fail;
428 #else /* insert virtual resistance */
429  VIRTUAL_RES ("no pivot != 0 found during Doolittle LU decomposition", c);
430 #endif
431  }
432 
433  // swap matrix rows if necessary and remember that step in the
434  // exchange table
435  if (c != pivot) {
436  A->exchangeRows (c, pivot);
437  Swap (int, rMap[c], rMap[pivot]);
438  Swap (nr_double_t, nPvt[c], nPvt[pivot]);
439  }
440 
441  // finally divide by the pivot element
442  if (c < N - 1) {
443  f = 1.0 / A_(c, c);
444  for (r = c + 1; r < N; r++) A_(r, c) *= f;
445  }
446  }
447 #if LU_FAILURE
448  fail:
449 #endif
450 }
451 
455 template <class nr_type_t>
457  nr_type_t f;
458  int i, c;
459 
460  // forward substitution in order to solve LY = B
461  for (i = 0; i < N; i++) {
462  f = B_(rMap[i]);
463  for (c = 0; c < i; c++) f -= A_(i, c) * X_(c);
464  X_(i) = f / A_(i, i);
465  }
466 
467  // backward substitution in order to solve UX = Y
468  for (i = N - 1; i >= 0; i--) {
469  f = X_(i);
470  for (c = i + 1; c < N; c++) f -= A_(i, c) * X_(c);
471  // remember that the Uii diagonal are ones only in Crout's definition
472  X_(i) = f;
473  }
474 }
475 
480 template <class nr_type_t>
482  nr_type_t f;
483  int i, c;
484 
485  // forward substitution in order to solve LY = B
486  for (i = 0; i < N; i++) {
487  f = B_(rMap[i]);
488  for (c = 0; c < i; c++) f -= A_(i, c) * X_(c);
489  // remember that the Lii diagonal are ones only in Doolittle's definition
490  X_(i) = f;
491  }
492 
493  // backward substitution in order to solve UX = Y
494  for (i = N - 1; i >= 0; i--) {
495  f = X_(i);
496  for (c = i + 1; c < N; c++) f -= A_(i, c) * X_(c);
497  X_(i) = f / A_(i, i);
498  }
499 }
500 
508 template <class nr_type_t>
510  nr_type_t f;
511  int error, conv, i, c, r;
512  int MaxIter = N; // -> less than N^3 operations
513  nr_double_t reltol = 1e-4;
514  nr_double_t abstol = NR_TINY;
515  nr_double_t diff, crit;
516 
517  // ensure that all diagonal values are non-zero
518  ensure_diagonal ();
519 
520  // try to raise diagonal dominance
521  preconditioner ();
522 
523  // decide here about possible convergence
524  if ((crit = convergence_criteria ()) >= 1) {
525 #if DEBUG && 0
526  logprint (LOG_STATUS, "NOTIFY: convergence criteria: %g >= 1 (%dx%d)\n",
527  crit, N, N);
528 #endif
529  //solve_lu ();
530  //return;
531  }
532 
533  // normalize the equation system to have ones on its diagonal
534  for (r = 0; r < N; r++) {
535  f = A_(r, r);
536  assert (f != 0); // singular matrix
537  for (c = 0; c < N; c++) A_(r, c) /= f;
538  B_(r) /= f;
539  }
540 
541  // the current X vector is a good initial guess for the iteration
542  tvector<nr_type_t> * Xprev = new tvector<nr_type_t> (*X);
543 
544  // start iterating here
545  i = 0; error = 0;
546  do {
547  // compute new solution vector
548  for (r = 0; r < N; r++) {
549  for (f = 0, c = 0; c < N; c++) {
550  if (algo == ALGO_GAUSS_SEIDEL) {
551  // Gauss-Seidel
552  if (c < r) f += A_(r, c) * X_(c);
553  else if (c > r) f += A_(r, c) * Xprev->get (c);
554  }
555  else {
556  // Jacobi
557  if (c != r) f += A_(r, c) * Xprev->get (c);
558  }
559  }
560  X_(r) = B_(r) - f;
561  }
562  // check for convergence
563  for (conv = 1, r = 0; r < N; r++) {
564  diff = abs (X_(r) - Xprev->get (r));
565  if (diff >= abstol + reltol * abs (X_(r))) {
566  conv = 0;
567  break;
568  }
569  if (!std::isfinite (diff)) { error++; break; }
570  }
571  // save last values
572  *Xprev = *X;
573  }
574  while (++i < MaxIter && !conv);
575 
576  delete Xprev;
577 
578  if (!conv || error) {
580  "WARNING: no convergence after %d %s iterations\n",
581  i, algo == ALGO_JACOBI ? "jacobi" : "gauss-seidel");
582  solve_lu_crout ();
583  }
584 #if DEBUG && 0
585  else {
587  "NOTIFY: %s convergence after %d iterations\n",
588  algo == ALGO_JACOBI ? "jacobi" : "gauss-seidel", i);
589  }
590 #endif
591 }
592 
597 template <class nr_type_t>
599  nr_type_t f;
600  int error, conv, i, c, r;
601  int MaxIter = N; // -> less than N^3 operations
602  nr_double_t reltol = 1e-4;
603  nr_double_t abstol = NR_TINY;
604  nr_double_t diff, crit, l = 1, d, s;
605 
606  // ensure that all diagonal values are non-zero
607  ensure_diagonal ();
608 
609  // try to raise diagonal dominance
610  preconditioner ();
611 
612  // decide here about possible convergence
613  if ((crit = convergence_criteria ()) >= 1) {
614 #if DEBUG && 0
615  logprint (LOG_STATUS, "NOTIFY: convergence criteria: %g >= 1 (%dx%d)\n",
616  crit, N, N);
617 #endif
618  //solve_lu ();
619  //return;
620  }
621 
622  // normalize the equation system to have ones on its diagonal
623  for (r = 0; r < N; r++) {
624  f = A_(r, r);
625  assert (f != 0); // singular matrix
626  for (c = 0; c < N; c++) A_(r, c) /= f;
627  B_(r) /= f;
628  }
629 
630  // the current X vector is a good initial guess for the iteration
631  tvector<nr_type_t> * Xprev = new tvector<nr_type_t> (*X);
632 
633  // start iterating here
634  i = 0; error = 0;
635  do {
636  // compute new solution vector
637  for (r = 0; r < N; r++) {
638  for (f = 0, c = 0; c < N; c++) {
639  if (c < r) f += A_(r, c) * X_(c);
640  else if (c > r) f += A_(r, c) * Xprev->get (c);
641  }
642  X_(r) = (1 - l) * Xprev->get (r) + l * (B_(r) - f);
643  }
644  // check for convergence
645  for (s = 0, d = 0, conv = 1, r = 0; r < N; r++) {
646  diff = abs (X_(r) - Xprev->get (r));
647  if (diff >= abstol + reltol * abs (X_(r))) {
648  conv = 0;
649  break;
650  }
651  d += diff; s += abs (X_(r));
652  if (!std::isfinite (diff)) { error++; break; }
653  }
654  if (!error) {
655  // adjust relaxation based on average errors
656  if ((s == 0 && d == 0) || d >= abstol * N + reltol * s) {
657  // values <= 1 -> non-convergence to convergence
658  if (l >= 0.6) l -= 0.1;
659  if (l >= 1.0) l = 1.0;
660  }
661  else {
662  // values >= 1 -> faster convergence
663  if (l < 1.5) l += 0.01;
664  if (l < 1.0) l = 1.0;
665  }
666  }
667  // save last values
668  *Xprev = *X;
669  }
670  while (++i < MaxIter && !conv);
671 
672  delete Xprev;
673 
674  if (!conv || error) {
676  "WARNING: no convergence after %d sor iterations (l = %g)\n",
677  i, l);
678  solve_lu_crout ();
679  }
680 #if DEBUG && 0
681  else {
683  "NOTIFY: sor convergence after %d iterations\n", i);
684  }
685 #endif
686 }
687 
691 template <class nr_type_t>
693  nr_double_t f = 0;
694  for (int r = 0; r < A->getCols (); r++) {
695  for (int c = 0; c < A->getCols (); c++) {
696  if (r != c) f += norm (A_(r, c) / A_(r, r));
697  }
698  }
699  return sqrt (f);
700 }
701 
705 template <class nr_type_t>
707  ensure_diagonal_MNA ();
708 }
709 
716 template <class nr_type_t>
718  int restart, exchanged, begin = 0, pairs;
719  int pivot1, pivot2, i;
720  do {
721  restart = exchanged = 0;
722  /* search for zero diagonals with lone pairs */
723  for (i = begin; i < N; i++) {
724  if (A_(i, i) == 0) {
725  pairs = countPairs (i, pivot1, pivot2);
726  if (pairs == 1) { /* lone pair found, substitute rows */
727  A->exchangeRows (i, pivot1);
728  B->exchangeRows (i, pivot1);
729  exchanged = 1;
730  }
731  else if ((pairs > 1) && !restart) {
732  restart = 1;
733  begin = i;
734  }
735  }
736  }
737 
738  /* all lone pairs are gone, look for zero diagonals with multiple pairs */
739  if (restart) {
740  for (i = begin; !exchanged && i < N; i++) {
741  if (A_(i, i) == 0) {
742  pairs = countPairs (i, pivot1, pivot2);
743  A->exchangeRows (i, pivot1);
744  B->exchangeRows (i, pivot1);
745  exchanged = 1;
746  }
747  }
748  }
749  }
750  while (restart);
751 }
752 
755 template <class nr_type_t>
756 int eqnsys<nr_type_t>::countPairs (int i, int& r1, int& r2) {
757  int pairs = 0;
758  for (int r = 0; r < N; r++) {
759  if (fabs (real (A_(r, i))) == 1.0) {
760  r1 = r;
761  pairs++;
762  for (r++; r < N; r++) {
763  if (fabs (real (A_(r, i))) == 1.0) {
764  r2 = r;
765  if (++pairs >= 2) return pairs;
766  }
767  }
768  }
769  }
770  return pairs;
771 }
772 
776 template <class nr_type_t>
778  int pivot, r;
779  nr_double_t MaxPivot;
780  for (int i = 0; i < N; i++) {
781  // find maximum column value for pivoting
782  for (MaxPivot = 0, pivot = i, r = 0; r < N; r++) {
783  if (abs (A_(r, i)) > MaxPivot &&
784  abs (A_(i, r)) >= abs (A_(r, r))) {
785  MaxPivot = abs (A_(r, i));
786  pivot = r;
787  }
788  }
789  // swap matrix rows if possible
790  if (i != pivot) {
791  A->exchangeRows (i, pivot);
792  B->exchangeRows (i, pivot);
793  }
794  }
795 }
796 
797 #define R_ (*R)
798 #define T_ (*T)
799 
802 template <class nr_type_t>
804  factorize_qrh ();
805  substitute_qrh ();
806 }
807 
810 template <class nr_type_t>
812  factorize_qr_householder ();
813  substitute_qr_householder ();
814 }
815 
819 template <class nr_type_t>
821  A->transpose ();
822  factorize_qr_householder ();
823  substitute_qr_householder_ls ();
824 }
825 
827 static void
828 euclidian_update (nr_double_t a, nr_double_t& n, nr_double_t& scale) {
829  nr_double_t x, ax;
830  if ((x = a) != 0) {
831  ax = fabs (x);
832  if (scale < ax) {
833  x = scale / ax;
834  n = 1 + n * x * x;
835  scale = ax;
836  }
837  else {
838  x = ax / scale;
839  n += x * x;
840  }
841  }
842 }
843 
846 template <class nr_type_t>
847 nr_double_t eqnsys<nr_type_t>::euclidian_c (int c, int r) {
848  nr_double_t scale = 0, n = 1;
849  for (int i = r; i < N; i++) {
850  euclidian_update (real (A_(i, c)), n, scale);
851  euclidian_update (imag (A_(i, c)), n, scale);
852  }
853  return scale * sqrt (n);
854 }
855 
858 template <class nr_type_t>
859 nr_double_t eqnsys<nr_type_t>::euclidian_r (int r, int c) {
860  nr_double_t scale = 0, n = 1;
861  for (int i = c; i < N; i++) {
862  euclidian_update (real (A_(r, i)), n, scale);
863  euclidian_update (imag (A_(r, i)), n, scale);
864  }
865  return scale * sqrt (n);
866 }
867 
868 template <typename nr_type_t>
869 inline nr_type_t cond_conj (nr_type_t t) {
870  return std::conj(t);
871 }
872 
873 template <>
874 inline double cond_conj (double t) {
875  return t;
876 }
877 
883 template <class nr_type_t>
885  int c, r, k, pivot;
886  nr_type_t f, t;
887  nr_double_t s, MaxPivot;
888 
889  if (R) delete R; R = new tvector<nr_type_t> (N);
890 
891  for (c = 0; c < N; c++) {
892  // compute column norms and save in work array
893  nPvt[c] = euclidian_c (c);
894  cMap[c] = c; // initialize permutation vector
895  }
896 
897  for (c = 0; c < N; c++) {
898 
899  // put column with largest norm into pivot position
900  MaxPivot = nPvt[c]; pivot = c;
901  for (r = c + 1; r < N; r++) {
902  if ((s = nPvt[r]) > MaxPivot) {
903  pivot = r;
904  MaxPivot = s;
905  }
906  }
907  if (pivot != c) {
908  A->exchangeCols (pivot, c);
909  Swap (int, cMap[pivot], cMap[c]);
910  Swap (nr_double_t, nPvt[pivot], nPvt[c]);
911  }
912 
913  // compute householder vector
914  if (c < N) {
915  nr_type_t a, b;
916  s = euclidian_c (c, c + 1);
917  a = A_(c, c);
918  b = -sign (a) * xhypot (a, s); // Wj
919  t = xhypot (s, a - b); // || Vi - Wi ||
920  R_(c) = b;
921  // householder vector entries Ui
922  A_(c, c) = (a - b) / t;
923  for (r = c + 1; r < N; r++) A_(r, c) /= t;
924  }
925  else {
926  R_(c) = A_(c, c);
927  }
928 
929  // apply householder transformation to remaining columns
930  for (r = c + 1; r < N; r++) {
931  for (f = 0, k = c; k < N; k++) f += cond_conj (A_(k, c)) * A_(k, r);
932  for (k = c; k < N; k++) A_(k, r) -= 2.0 * f * A_(k, c);
933  }
934 
935  // update norms of remaining columns too
936  for (r = c + 1; r < N; r++) {
937  nPvt[r] = euclidian_c (r, c + 1);
938  }
939  }
940 }
941 
948 template <class nr_type_t>
950  int c, r, pivot;
951  nr_double_t s, MaxPivot;
952 
953  if (T) delete T; T = new tvector<nr_type_t> (N);
954 
955  for (c = 0; c < N; c++) {
956  // compute column norms and save in work array
957  nPvt[c] = euclidian_c (c);
958  cMap[c] = c; // initialize permutation vector
959  }
960 
961  for (c = 0; c < N; c++) {
962 
963  // put column with largest norm into pivot position
964  MaxPivot = nPvt[c]; pivot = c;
965  for (r = c + 1; r < N; r++)
966  if ((s = nPvt[r]) > MaxPivot) {
967  pivot = r;
968  MaxPivot = s;
969  }
970  if (pivot != c) {
971  A->exchangeCols (pivot, c);
972  Swap (int, cMap[pivot], cMap[c]);
973  Swap (nr_double_t, nPvt[pivot], nPvt[c]);
974  }
975 
976  // compute and apply householder vector
977  T_(c) = householder_left (c);
978 
979  // update norms of remaining columns too
980  for (r = c + 1; r < N; r++) {
981  if ((s = nPvt[r]) > 0) {
982  nr_double_t y = 0;
983  nr_double_t t = norm (A_(c, r) / s);
984  if (t < 1)
985  y = s * sqrt (1 - t);
986  if (fabs (y / s) < NR_TINY)
987  nPvt[r] = euclidian_c (r, c + 1);
988  else
989  nPvt[r] = y;
990  }
991  }
992  }
993 }
994 
997 template <class nr_type_t>
999  int c, r;
1000  nr_type_t f;
1001 
1002  // form the new right hand side Q'B
1003  for (c = 0; c < N - 1; c++) {
1004  // scalar product u_k^T * B
1005  for (f = 0, r = c; r < N; r++) f += cond_conj (A_(r, c)) * B_(r);
1006  // z - 2 * f * u_k
1007  for (r = c; r < N; r++) B_(r) -= 2.0 * f * A_(r, c);
1008  }
1009 
1010  // backward substitution in order to solve RX = Q'B
1011  for (r = N - 1; r >= 0; r--) {
1012  f = B_(r);
1013  for (c = r + 1; c < N; c++) f -= A_(r, c) * X_(cMap[c]);
1014  if (abs (R_(r)) > NR_EPSI)
1015  X_(cMap[r]) = f / R_(r);
1016  else
1017  X_(cMap[r]) = 0;
1018  }
1019 }
1020 
1023 template <class nr_type_t>
1025  int c, r;
1026  nr_type_t f;
1027 
1028  // form the new right hand side Q'B
1029  for (c = 0; c < N; c++) {
1030  if (T_(c) != 0) {
1031  // scalar product u' * B
1032  for (f = B_(c), r = c + 1; r < N; r++) f += cond_conj (A_(r, c)) * B_(r);
1033  // z - T * f * u
1034  f *= cond_conj (T_(c)); B_(c) -= f;
1035  for (r = c + 1; r < N; r++) B_(r) -= f * A_(r, c);
1036  }
1037  }
1038 
1039  // backward substitution in order to solve RX = Q'B
1040  for (r = N - 1; r >= 0; r--) {
1041  for (f = B_(r), c = r + 1; c < N; c++) f -= A_(r, c) * X_(cMap[c]);
1042  if (abs (A_(r, r)) > NR_EPSI)
1043  X_(cMap[r]) = f / A_(r, r);
1044  else
1045  X_(cMap[r]) = 0;
1046  }
1047 }
1048 
1053 template <class nr_type_t>
1055  int c, r;
1056  nr_type_t f;
1057 
1058  // forward substitution in order to solve R'X = B
1059  for (r = 0; r < N; r++) {
1060  for (f = B_(r), c = 0; c < r; c++) f -= A_(c, r) * B_(c);
1061  if (abs (A_(r, r)) > NR_EPSI)
1062  B_(r) = f / A_(r, r);
1063  else
1064  B_(r) = 0;
1065  }
1066 
1067  // compute the least square solution QX
1068  for (c = N - 1; c >= 0; c--) {
1069  if (T_(c) != 0) {
1070  // scalar product u' * B
1071  for (f = B_(c), r = c + 1; r < N; r++) f += cond_conj (A_(r, c)) * B_(r);
1072  // z - T * f * u_k
1073  f *= T_(c); B_(c) -= f;
1074  for (r = c + 1; r < N; r++) B_(r) -= f * A_(r, c);
1075  }
1076  }
1077 
1078  // permute solution vector
1079  for (r = 0; r < N; r++) X_(cMap[r]) = B_(r);
1080 }
1081 
1082 #define sign_(a) (real (a) < 0 ? -1 : 1)
1083 
1090 template <class nr_type_t>
1092  nr_type_t a, b, t;
1093  nr_double_t s, g;
1094  s = euclidian_c (c, c + 1);
1095  if (s == 0 && imag (A_(c, c)) == 0) {
1096  // no reflection necessary
1097  t = 0;
1098  }
1099  else {
1100  // calculate householder reflection
1101  a = A_(c, c);
1102  g = sign_(a) * xhypot (a, s);
1103  b = a + g;
1104  t = b / g;
1105  // store householder vector
1106  for (int r = c + 1; r < N; r++) A_(r, c) /= b;
1107  A_(c, c) = -g;
1108  }
1109  return t;
1110 }
1111 
1116 template <class nr_type_t>
1118  // compute householder vector
1119  nr_type_t t = householder_create_left (c);
1120  // apply householder transformation to remaining columns if necessary
1121  if (t != 0) {
1122  householder_apply_left (c, t);
1123  }
1124  return t;
1125 }
1126 
1131 template <class nr_type_t>
1133  // compute householder vector
1134  nr_type_t t = householder_create_right (r);
1135  // apply householder transformation to remaining rows if necessary
1136  if (t != 0) {
1137  householder_apply_right (r, t);
1138  }
1139  return t;
1140 }
1141 
1148 template <class nr_type_t>
1150  nr_type_t a, b, t;
1151  nr_double_t s, g;
1152  s = euclidian_r (r, r + 2);
1153  if (s == 0 && imag (A_(r, r + 1)) == 0) {
1154  // no reflection necessary
1155  t = 0;
1156  }
1157  else {
1158  // calculate householder reflection
1159  a = A_(r, r + 1);
1160  g = sign_(a) * xhypot (a, s);
1161  b = a + g;
1162  t = b / g;
1163  // store householder vector
1164  for (int c = r + 2; c < N; c++) A_(r, c) /= b;
1165  A_(r, r + 1) = -g;
1166  }
1167  return t;
1168 }
1169 
1172 template <class nr_type_t>
1174  nr_type_t f;
1175  int k, r;
1176  // apply the householder vector to each right-hand column
1177  for (r = c + 1; r < N; r++) {
1178  // calculate f = u' * A (a scalar product)
1179  f = A_(c, r);
1180  for (k = c + 1; k < N; k++) f += cond_conj (A_(k, c)) * A_(k, r);
1181  // calculate A -= T * f * u
1182  f *= cond_conj (t); A_(c, r) -= f;
1183  for (k = c + 1; k < N; k++) A_(k, r) -= f * A_(k, c);
1184  }
1185 }
1186 
1189 template <class nr_type_t>
1191  nr_type_t f;
1192  int k, c;
1193  // apply the householder vector to each downward row
1194  for (c = r + 1; c < N; c++) {
1195  // calculate f = u' * A (a scalar product)
1196  f = A_(c, r + 1);
1197  for (k = r + 2; k < N; k++) f += cond_conj (A_(r, k)) * A_(c, k);
1198  // calculate A -= T * f * u
1199  f *= cond_conj (t); A_(c, r + 1) -= f;
1200  for (k = r + 2; k < N; k++) A_(c, k) -= f * A_(r, k);
1201  }
1202 }
1203 
1204 // Some helper defines for SVD.
1205 #define V_ (*V)
1206 #define S_ (*S)
1207 #define E_ (*E)
1208 #define U_ (*A)
1209 
1213 template <class nr_type_t>
1215  nr_type_t f;
1216  int k, c;
1217  // apply the householder vector to each downward row
1218  for (c = r + 1; c < N; c++) {
1219  // calculate f = u' * A (a scalar product)
1220  f = V_(c, r + 1);
1221  for (k = r + 2; k < N; k++) f += cond_conj (A_(r, k)) * V_(c, k);
1222  // calculate A -= T * f * u
1223  f *= cond_conj (t); V_(c, r + 1) -= f;
1224  for (k = r + 2; k < N; k++) V_(c, k) -= f * A_(r, k);
1225  }
1226 }
1227 
1230 template <class nr_type_t>
1232  factorize_svd ();
1233  chop_svd ();
1234  substitute_svd ();
1235 }
1236 
1238 template <class nr_type_t>
1240  int c;
1241  nr_double_t Max, Min;
1242  Max = 0.0;
1243  for (c = 0; c < N; c++) if (fabs (S_(c)) > Max) Max = fabs (S_(c));
1244  Min = Max * NR_EPSI;
1245  for (c = 0; c < N; c++) if (fabs (S_(c)) < Min) S_(c) = 0.0;
1246 }
1247 
1251 template <class nr_type_t>
1253  int c, r;
1254  nr_type_t f;
1255  // calculate U'B
1256  for (c = 0; c < N; c++) {
1257  f = 0.0;
1258  // non-zero result only if S is non-zero
1259  if (S_(c) != 0.0) {
1260  for (r = 0; r < N; r++) f += cond_conj (U_(r, c)) * B_(r);
1261  // this is the divide by S
1262  f /= S_(c);
1263  }
1264  R_(c) = f;
1265  }
1266  // matrix multiply by V to get the final solution
1267  for (r = 0; r < N; r++) {
1268  for (f = 0.0, c = 0; c < N; c++) f += cond_conj (V_(c, r)) * R_(c);
1269  X_(r) = f;
1270  }
1271 }
1272 
1277 template <class nr_type_t>
1279  int i, j, l;
1280  nr_type_t t;
1281 
1282  // allocate space for vectors and matrices
1283  if (R) delete R; R = new tvector<nr_type_t> (N);
1284  if (T) delete T; T = new tvector<nr_type_t> (N);
1285  if (V) delete V; V = new tmatrix<nr_type_t> (N);
1286  if (S) delete S; S = new tvector<nr_double_t> (N);
1287  if (E) delete E; E = new tvector<nr_double_t> (N);
1288 
1289  // bidiagonalization through householder transformations
1290  for (i = 0; i < N; i++) {
1291  T_(i) = householder_left (i);
1292  if (i < N - 1) R_(i) = householder_right (i);
1293  }
1294 
1295  // copy over the real valued bidiagonal values
1296  for (i = 0; i < N; i++) S_(i) = real (A_(i, i));
1297  for (E_(0) = 0, i = 1; i < N; i++) E_(i) = real (A_(i - 1, i));
1298 
1299  // backward accumulation of right-hand householder transformations
1300  // yields the V' matrix
1301  for (l = N, i = N - 1; i >= 0; l = i--) {
1302  if (i < N - 1) {
1303  if ((t = R_(i)) != 0.0) {
1304  householder_apply_right_extern (i, cond_conj (t));
1305  }
1306  else for (j = l; j < N; j++) // cleanup this row
1307  V_(i, j) = V_(j, i) = 0.0;
1308  }
1309  V_(i, i) = 1.0;
1310  }
1311 
1312  // backward accumulation of left-hand householder transformations
1313  // yields the U matrix in place of the A matrix
1314  for (l = N, i = N - 1; i >= 0; l = i--) {
1315  for (j = l; j < N; j++) // cleanup upper row
1316  A_(i, j) = 0.0;
1317  if ((t = T_(i)) != 0.0) {
1318  householder_apply_left (i, cond_conj (t));
1319  for (j = l; j < N; j++) A_(j, i) *= -t;
1320  }
1321  else for (j = l; j < N; j++) // cleanup this column
1322  A_(j, i) = 0.0;
1323  A_(i, i) = 1.0 - t;
1324  }
1325 
1326  // S and E contain diagonal and super-diagonal, A contains U, V'
1327  // calculated; now diagonalization can begin
1328  diagonalize_svd ();
1329 }
1330 
1331 #ifndef MAX
1332 # define MAX(x,y) (((x) > (y)) ? (x) : (y))
1333 #endif
1334 
1336 static nr_double_t
1337 givens (nr_double_t a, nr_double_t b, nr_double_t& c, nr_double_t& s) {
1338  nr_double_t z = xhypot (a, b);
1339  c = a / z;
1340  s = b / z;
1341  return z;
1342 }
1343 
1344 template <class nr_type_t>
1346  nr_double_t c, nr_double_t s) {
1347  for (int i = 0; i < N; i++) {
1348  nr_type_t y = U_(i, c1);
1349  nr_type_t z = U_(i, c2);
1350  U_(i, c1) = y * c + z * s;
1351  U_(i, c2) = z * c - y * s;
1352  }
1353 }
1354 
1355 template <class nr_type_t>
1357  nr_double_t c, nr_double_t s) {
1358  for (int i = 0; i < N; i++) {
1359  nr_type_t y = V_(r1, i);
1360  nr_type_t z = V_(r2, i);
1361  V_(r1, i) = y * c + z * s;
1362  V_(r2, i) = z * c - y * s;
1363  }
1364 }
1365 
1371 template <class nr_type_t>
1373  bool split;
1374  int i, l, j, its, k, n, MaxIters = 30;
1375  nr_double_t an, f, g, h, d, c, s, b, a;
1376 
1377  // find largest bidiagonal value
1378  for (an = 0, i = 0; i < N; i++)
1379  an = MAX (an, fabs (S_(i)) + fabs (E_(i)));
1380 
1381  // diagonalize the bidiagonal matrix (stored as super-diagonal
1382  // vector E and diagonal vector S)
1383  for (k = N - 1; k >= 0; k--) {
1384  // loop over singular values
1385  for (its = 0; its <= MaxIters; its++) {
1386  split = true;
1387  // check for a zero entry along the super-diagonal E, if there
1388  // is one, it is possible to QR iterate on two separate matrices
1389  // above and below it
1390  for (n = 0, l = k; l >= 1; l--) {
1391  // note that E_(0) is always zero
1392  n = l - 1;
1393  if (fabs (E_(l)) + an == an) { split = false; break; }
1394  if (fabs (S_(n)) + an == an) break;
1395  }
1396  // if there is a zero on the diagonal S, it is possible to zero
1397  // out the corresponding super-diagonal E entry to its right
1398  if (split) {
1399  // cancellation of E_(l), if l > 0
1400  c = 0.0;
1401  s = 1.0;
1402  for (i = l; i <= k; i++) {
1403  f = -s * E_(i);
1404  E_(i) *= c;
1405  if (fabs (f) + an == an) break;
1406  g = S_(i);
1407  S_(i) = givens (f, g, c, s);
1408  // apply inverse rotation to U
1409  givens_apply_u (n, i, c, s);
1410  }
1411  }
1412 
1413  d = S_(k);
1414  // convergence
1415  if (l == k) {
1416  // singular value is made non-negative
1417  if (d < 0.0) {
1418  S_(k) = -d;
1419  for (j = 0; j < N; j++) V_(k, j) = -V_(k, j);
1420  }
1421  break;
1422  }
1423  if (its == MaxIters) {
1424  logprint (LOG_ERROR, "WARNING: no convergence in %d SVD iterations\n",
1425  MaxIters);
1426  }
1427 
1428  // shift from bottom 2-by-2 minor
1429  a = S_(l);
1430  n = k - 1;
1431  b = S_(n);
1432  g = E_(n);
1433  h = E_(k);
1434 
1435  // compute QR shift value (as close as possible to the largest
1436  // eigenvalue of the 2-by-2 minor matrix)
1437  f = (b - d) * (b + d) + (g - h) * (g + h);
1438  f /= 2.0 * h * b;
1439  f += sign_(f) * xhypot (f, 1.0);
1440  f = ((a - d) * (a + d) + h * (b / f - h)) / a;
1441  // f => (B_{ll}^2 - u) / B_{ll}
1442  // u => eigenvalue of T = B' * B nearer T_{22} (Wilkinson shift)
1443 
1444  // next QR transformation
1445  c = s = 1.0;
1446  for (j = l; j <= n; j++) {
1447  i = j + 1;
1448  g = E_(i);
1449  b = S_(i);
1450  h = s * g; // h => right-hand non-zero to annihilate
1451  g *= c;
1452  E_(j) = givens (f, h, c, s);
1453  // perform the rotation
1454  f = a * c + g * s;
1455  g = g * c - a * s;
1456  h = b * s;
1457  b *= c;
1458  // here: +- -+
1459  // | f g | = B * V'_j (also first V'_1)
1460  // | h b |
1461  // +- -+
1462 
1463  // accumulate the rotation in V'
1464  givens_apply_v (j, i, c, s);
1465  d = S_(j) = xhypot (f, h);
1466  // rotation can be arbitrary if d = 0
1467  if (d != 0.0) {
1468  // d => non-zero result on diagonal
1469  d = 1.0 / d;
1470  // rotation coefficients to annihilate the lower non-zero
1471  c = f * d;
1472  s = h * d;
1473  }
1474  f = c * g + s * b;
1475  a = c * b - s * g;
1476  // here: +- -+
1477  // | d f | => U_j * B
1478  // | 0 a |
1479  // +- -+
1480 
1481  // accumulate rotation into U
1482  givens_apply_u (j, i, c, s);
1483  }
1484  E_(l) = 0;
1485  E_(k) = f;
1486  S_(k) = a;
1487  }
1488  }
1489 }
1490 
1491 } // namespace qucs
1492 
1493 #undef V_
matrix inverse(matrix a)
Compute inverse matrix.
Definition: matrix.cpp:847
l
Definition: parse_vcd.y:213
#define R_
Definition: eqnsys.cpp:797
void ensure_diagonal_MNA(void)
Definition: eqnsys.cpp:717
#define S_
Definition: eqnsys.cpp:1206
void solve_inverse(void)
Definition: eqnsys.cpp:183
void setData(int d)
Definition: exception.h:56
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 N(n)
Definition: equation.cpp:58
tvector< nr_type_t > * X
Definition: eqnsys.h:87
void householder_apply_left(int, nr_type_t)
Definition: eqnsys.cpp:1173
#define sign_(a)
Definition: eqnsys.cpp:1082
void solve_gauss_jordan(void)
Definition: eqnsys.cpp:234
void ensure_diagonal(void)
Definition: eqnsys.cpp:706
void solve_lu_crout(void)
Definition: eqnsys.cpp:288
eqnsys()
Constructor creates an unnamed instance of the eqnsys class.
Definition: eqnsys.cpp:51
nr_double_t convergence_criteria(void)
Definition: eqnsys.cpp:692
nr_type_t householder_left(int)
Definition: eqnsys.cpp:1117
void factorize_lu_crout(void)
Definition: eqnsys.cpp:319
void factorize_qrh(void)
Definition: eqnsys.cpp:884
nr_type_t householder_create_left(int)
Definition: eqnsys.cpp:1091
t
Definition: parse_vcd.y:290
void substitute_lu_doolittle(void)
Definition: eqnsys.cpp:481
#define E(equ)
Definition: parasweep.cpp:81
#define V_
Definition: eqnsys.cpp:1205
void substitute_qr_householder_ls(void)
Definition: eqnsys.cpp:1054
static nr_double_t givens(nr_double_t a, nr_double_t b, nr_double_t &c, nr_double_t &s)
Helper function computes Givens rotation.
Definition: eqnsys.cpp:1337
void passEquationSys(tmatrix< nr_type_t > *, tvector< nr_type_t > *, tvector< nr_type_t > *)
Definition: eqnsys.cpp:99
void householder_apply_right(int, nr_type_t)
Definition: eqnsys.cpp:1190
void factorize_qr_householder(void)
Definition: eqnsys.cpp:949
n
Definition: parse_citi.y:147
#define S(con)
Definition: property.cpp:158
void factorize_lu_doolittle(void)
Definition: eqnsys.cpp:385
h
Definition: parse_vcd.y:214
r
Definition: parse_mdl.y:515
nr_complex_t sign(const nr_complex_t z)
complex sign function
Definition: complex.cpp:435
#define throw_exception(e)
#define R(con)
i
Definition: parse_mdl.y:516
matrix imag(matrix a)
Imaginary part matrix.
Definition: matrix.cpp:581
void substitute_qrh(void)
Definition: eqnsys.cpp:998
nr_complex_t sqrt(const nr_complex_t z)
Compute principal value of square root.
Definition: complex.cpp:271
nr_double_t xhypot(const nr_complex_t a, const nr_complex_t b)
Euclidean distance function for complex argument.
Definition: complex.cpp:465
#define T_
Definition: eqnsys.cpp:798
void solve_svd(void)
Definition: eqnsys.cpp:1231
void givens_apply_v(int, int, nr_double_t, nr_double_t)
Definition: eqnsys.cpp:1356
#define B_
Definition: eqnsys.cpp:189
void substitute_svd(void)
Definition: eqnsys.cpp:1252
void solve_lu_doolittle(void)
Definition: eqnsys.cpp:302
nr_type_t get(int)
Definition: tvector.cpp:101
void factorize_svd(void)
Definition: eqnsys.cpp:1278
#define V(con)
Definition: evaluate.cpp:63
eqn::constant * c1
#define X_
Definition: eqnsys.cpp:188
~eqnsys()
Destructor deletes the eqnsys class object.
Definition: eqnsys.cpp:65
eqn::constant * c2
#define A_
Definition: eqnsys.cpp:187
x
Definition: parse_mdl.y:498
static void euclidian_update(nr_double_t a, nr_double_t &n, nr_double_t &scale)
Definition: eqnsys.cpp:828
#define Swap(type, a, b)
Little helper macro.
Definition: eqnsys.cpp:45
void preconditioner(void)
Definition: eqnsys.cpp:777
tvector< nr_type_t > * B
Definition: eqnsys.h:86
void solve_qr_ls(void)
Definition: eqnsys.cpp:820
void solve(void)
Definition: eqnsys.cpp:124
nr_double_t euclidian_r(int, int c=1)
Definition: eqnsys.cpp:859
void chop_svd(void)
Annihilates near-zero singular values.
Definition: eqnsys.cpp:1239
void solve_gauss(void)
Definition: eqnsys.cpp:194
vector diff(vector, vector, int n=1)
Definition: vector.cpp:591
void householder_apply_right_extern(int, nr_type_t)
Definition: eqnsys.cpp:1214
nr_type_t cond_conj(nr_type_t t)
Definition: eqnsys.cpp:869
void substitute_qr_householder(void)
Definition: eqnsys.cpp:1024
void solve_qr(void)
Definition: eqnsys.cpp:811
#define U_
Definition: eqnsys.cpp:1208
nr_double_t euclidian_c(int, int r=1)
Definition: eqnsys.cpp:847
#define B(con)
Definition: evaluate.cpp:70
nr_double_t norm(const nr_complex_t z)
Compute euclidian norm of complex number.
Definition: complex.cpp:283
void solve_sor(void)
Definition: eqnsys.cpp:598
pairs
tmatrix< nr_type_t > * A
Definition: eqnsys.h:84
#define A(a)
Definition: eqndefined.cpp:72
void diagonalize_svd(void)
Definition: eqnsys.cpp:1372
int countPairs(int, int &, int &)
Definition: eqnsys.cpp:756
y
Definition: parse_mdl.y:499
void substitute_lu_crout(void)
Definition: eqnsys.cpp:456
void setText(const char *,...)
Definition: exception.cpp:67
nr_type_t householder_create_right(int)
Definition: eqnsys.cpp:1149
#define LOG_ERROR
Definition: logging.h:28
matrix conj(matrix a)
Conjugate complex matrix.
Definition: matrix.cpp:505
#define LOG_STATUS
Definition: logging.h:29
#define VIRTUAL_RES(txt, i)
Definition: eqnsys.cpp:275
vcd scale
Definition: parse_vcd.y:180
void logprint(int level, const char *format,...)
Definition: logging.c:37
void solve_qrh(void)
Definition: eqnsys.cpp:803
void givens_apply_u(int, int, nr_double_t, nr_double_t)
Definition: eqnsys.cpp:1345
#define E_
Definition: eqnsys.cpp:1207
nr_type_t householder_right(int)
Definition: eqnsys.cpp:1132
void solve_iterative(void)
Definition: eqnsys.cpp:509
#define MAX(x, y)
Maximum of x and y.
Definition: constants.h:127