28 #include "qucs_typedefs.h"
45 #define Swap(type,a,b) { type t; t = a; a = b; b = t; }
50 template <
class nr_type_t>
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;
79 template <
class nr_type_t>
98 template <
class nr_type_t>
105 if (
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];
115 if (
B != NULL)
delete B;
123 template <
class nr_type_t>
126 time_t
t = time (NULL);
136 solve_gauss_jordan ();
142 solve_lu_doolittle ();
145 factorize_lu_crout ();
148 factorize_lu_doolittle ();
151 substitute_lu_crout ();
154 substitute_lu_doolittle ();
177 A->getRows (),
A->getCols (), time (NULL) -
t);
182 template <
class nr_type_t>
193 template <
class nr_type_t>
195 nr_double_t MaxPivot;
200 for (i = 0; i <
N; i++) {
202 for (MaxPivot = 0, pivot = r = i; r <
N; r++) {
203 if (
abs (
A_(r, i)) > MaxPivot) {
204 MaxPivot =
abs (
A_(r, i));
209 assert (MaxPivot != 0);
211 A->exchangeRows (i, pivot);
212 B->exchangeRows (i, pivot);
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);
223 for (i = N - 1; i >= 0; i--) {
225 for (c = i + 1; c <
N; c++) f -=
A_(i, c) *
X_(c);
226 X_(i) = f /
A_(i, i);
233 template <
class nr_type_t>
235 nr_double_t MaxPivot;
240 for (i = 0; i <
N; i++) {
242 for (MaxPivot = 0, pivot = r = i; r <
N; r++) {
243 if (
abs (
A_(r, i)) > MaxPivot) {
244 MaxPivot =
abs (
A_(r, i));
249 assert (MaxPivot != 0);
251 A->exchangeRows (i, pivot);
252 B->exchangeRows (i, pivot);
257 for (c = i + 1; c <
N; c++)
A_(i, c) /= f;
261 for (r = 0; r <
N; r++) {
264 for (c = i + 1; c <
N; c++)
A_(r, c) -= f *
A_(i, c);
275 #define VIRTUAL_RES(txt,i) { \
276 qucs::exception * e = new qucs::exception (EXCEPTION_SINGULAR); \
278 e->setData (rMap[i]); \
279 A_(i, i) = NR_TINY; \
280 throw_exception (e); }
287 template <
class nr_type_t>
293 factorize_lu_crout ();
297 substitute_lu_crout ();
301 template <
class nr_type_t>
307 factorize_lu_doolittle ();
311 substitute_lu_doolittle ();
318 template <
class nr_type_t>
320 nr_double_t
d, MaxPivot;
325 for (r = 0; r <
N; r++) {
326 for (MaxPivot = 0, c = 0; c <
N; c++)
327 if ((d =
abs (
A_(r, c))) > MaxPivot)
329 if (MaxPivot <= 0) MaxPivot = NR_TINY;
330 nPvt[
r] = 1 / MaxPivot;
335 for (c = 0; c <
N; c++) {
337 for (r = 0; r <
c; r++) {
339 for (k = 0; k <
r; k++) f -=
A_(r, k) *
A_(k, c);
340 A_(r, c) = f /
A_(r, r);
343 for (MaxPivot = 0, pivot = r; r <
N; r++) {
345 for (k = 0; k <
c; k++) f -=
A_(r, k) *
A_(k, c);
348 if ((d = nPvt[r] *
abs (f)) > MaxPivot) {
358 e->
setText (
"no pivot != 0 found during Crout LU decomposition");
363 VIRTUAL_RES (
"no pivot != 0 found during Crout LU decomposition", c);
370 A->exchangeRows (c, pivot);
371 Swap (
int, rMap[c], rMap[pivot]);
372 Swap (nr_double_t, nPvt[c], nPvt[pivot]);
384 template <
class nr_type_t>
386 nr_double_t
d, MaxPivot;
391 for (r = 0; r <
N; r++) {
392 for (MaxPivot = 0, c = 0; c <
N; c++)
393 if ((d =
abs (
A_(r, c))) > MaxPivot)
395 if (MaxPivot <= 0) MaxPivot = NR_TINY;
396 nPvt[
r] = 1 / MaxPivot;
401 for (c = 0; c <
N; c++) {
403 for (r = 0; r <
c; r++) {
405 for (k = 0; k <
r; k++) f -=
A_(r, k) *
A_(k, c);
409 for (MaxPivot = 0, pivot = r; r <
N; r++) {
411 for (k = 0; k <
c; k++) f -=
A_(r, k) *
A_(k, c);
414 if ((d = nPvt[r] *
abs (f)) > MaxPivot) {
424 e->
setText (
"no pivot != 0 found during Doolittle LU decomposition");
429 VIRTUAL_RES (
"no pivot != 0 found during Doolittle LU decomposition", c);
436 A->exchangeRows (c, pivot);
437 Swap (
int, rMap[c], rMap[pivot]);
438 Swap (nr_double_t, nPvt[c], nPvt[pivot]);
444 for (r = c + 1; r <
N; r++)
A_(r, c) *= f;
455 template <
class nr_type_t>
461 for (i = 0; i <
N; i++) {
463 for (c = 0; c <
i; c++) f -=
A_(i, c) *
X_(c);
464 X_(i) = f /
A_(i, i);
468 for (i = N - 1; i >= 0; i--) {
470 for (c = i + 1; c <
N; c++) f -=
A_(i, c) *
X_(c);
480 template <
class nr_type_t>
486 for (i = 0; i <
N; i++) {
488 for (c = 0; c <
i; c++) f -=
A_(i, c) *
X_(c);
494 for (i = N - 1; i >= 0; i--) {
496 for (c = i + 1; c <
N; c++) f -=
A_(i, c) *
X_(c);
497 X_(i) = f /
A_(i, i);
508 template <
class nr_type_t>
511 int error, conv,
i,
c,
r;
513 nr_double_t reltol = 1e-4;
514 nr_double_t abstol = NR_TINY;
515 nr_double_t
diff, crit;
524 if ((crit = convergence_criteria ()) >= 1) {
534 for (r = 0; r <
N; r++) {
537 for (c = 0; c <
N; c++)
A_(r, c) /= f;
548 for (r = 0; r <
N; r++) {
549 for (f = 0, c = 0; c <
N; c++) {
552 if (c < r) f +=
A_(r, c) *
X_(c);
553 else if (c > r) f +=
A_(r, c) * Xprev->
get (c);
557 if (c != r) f +=
A_(r, c) * Xprev->
get (c);
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))) {
569 if (!std::isfinite (diff)) { error++;
break; }
574 while (++i < MaxIter && !conv);
578 if (!conv || error) {
580 "WARNING: no convergence after %d %s iterations\n",
581 i, algo ==
ALGO_JACOBI ?
"jacobi" :
"gauss-seidel");
587 "NOTIFY: %s convergence after %d iterations\n",
588 algo ==
ALGO_JACOBI ?
"jacobi" :
"gauss-seidel", i);
597 template <
class nr_type_t>
600 int error, conv,
i,
c,
r;
602 nr_double_t reltol = 1e-4;
603 nr_double_t abstol = NR_TINY;
604 nr_double_t
diff, crit,
l = 1,
d,
s;
613 if ((crit = convergence_criteria ()) >= 1) {
623 for (r = 0; r <
N; r++) {
626 for (c = 0; c <
N; c++)
A_(r, c) /= f;
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);
642 X_(r) = (1 -
l) * Xprev->
get (r) + l * (
B_(r) - f);
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))) {
652 if (!std::isfinite (diff)) { error++;
break; }
656 if ((s == 0 &&
d == 0) ||
d >= abstol * N + reltol * s) {
658 if (l >= 0.6) l -= 0.1;
659 if (l >= 1.0) l = 1.0;
663 if (l < 1.5) l += 0.01;
664 if (l < 1.0) l = 1.0;
670 while (++i < MaxIter && !conv);
674 if (!conv || error) {
676 "WARNING: no convergence after %d sor iterations (l = %g)\n",
683 "NOTIFY: sor convergence after %d iterations\n", i);
691 template <
class nr_type_t>
694 for (
int r = 0;
r <
A->getCols ();
r++) {
695 for (
int c = 0;
c <
A->getCols ();
c++) {
705 template <
class nr_type_t>
707 ensure_diagonal_MNA ();
716 template <
class nr_type_t>
718 int restart, exchanged, begin = 0,
pairs;
719 int pivot1, pivot2,
i;
721 restart = exchanged = 0;
723 for (i = begin; i <
N; i++) {
725 pairs = countPairs (i, pivot1, pivot2);
727 A->exchangeRows (i, pivot1);
728 B->exchangeRows (i, pivot1);
731 else if ((
pairs > 1) && !restart) {
740 for (i = begin; !exchanged && i <
N; i++) {
742 pairs = countPairs (i, pivot1, pivot2);
743 A->exchangeRows (i, pivot1);
744 B->exchangeRows (i, pivot1);
755 template <
class nr_type_t>
758 for (
int r = 0;
r <
N;
r++) {
759 if (fabs (
real (
A_(
r, i))) == 1.0) {
762 for (
r++;
r <
N;
r++) {
763 if (fabs (
real (
A_(
r, i))) == 1.0) {
765 if (++pairs >= 2)
return pairs;
776 template <
class nr_type_t>
779 nr_double_t MaxPivot;
780 for (
int i = 0;
i <
N;
i++) {
782 for (MaxPivot = 0, pivot =
i, r = 0; r <
N; r++) {
783 if (
abs (
A_(r,
i)) > MaxPivot &&
785 MaxPivot =
abs (
A_(r,
i));
791 A->exchangeRows (
i, pivot);
792 B->exchangeRows (
i, pivot);
802 template <
class nr_type_t>
810 template <
class nr_type_t>
812 factorize_qr_householder ();
813 substitute_qr_householder ();
819 template <
class nr_type_t>
822 factorize_qr_householder ();
823 substitute_qr_householder_ls ();
846 template <
class nr_type_t>
848 nr_double_t
scale = 0,
n = 1;
849 for (
int i = r;
i <
N;
i++) {
853 return scale *
sqrt (
n);
858 template <
class nr_type_t>
860 nr_double_t
scale = 0,
n = 1;
861 for (
int i = c;
i <
N;
i++) {
865 return scale *
sqrt (
n);
868 template <
typename nr_type_t>
883 template <
class nr_type_t>
887 nr_double_t
s, MaxPivot;
891 for (c = 0; c <
N; c++) {
893 nPvt[
c] = euclidian_c (c);
897 for (c = 0; c <
N; c++) {
900 MaxPivot = nPvt[
c]; pivot =
c;
901 for (r = c + 1; r <
N; r++) {
902 if ((s = nPvt[r]) > MaxPivot) {
908 A->exchangeCols (pivot, c);
909 Swap (
int, cMap[pivot], cMap[c]);
910 Swap (nr_double_t, nPvt[pivot], nPvt[c]);
916 s = euclidian_c (c, c + 1);
922 A_(c, c) = (a - b) / t;
923 for (r = c + 1; r <
N; r++)
A_(r, c) /=
t;
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);
936 for (r = c + 1; r <
N; r++) {
937 nPvt[
r] = euclidian_c (r, c + 1);
948 template <
class nr_type_t>
951 nr_double_t
s, MaxPivot;
955 for (c = 0; c <
N; c++) {
957 nPvt[
c] = euclidian_c (c);
961 for (c = 0; c <
N; c++) {
964 MaxPivot = nPvt[
c]; pivot =
c;
965 for (r = c + 1; r <
N; r++)
966 if ((s = nPvt[r]) > MaxPivot) {
971 A->exchangeCols (pivot, c);
972 Swap (
int, cMap[pivot], cMap[c]);
973 Swap (nr_double_t, nPvt[pivot], nPvt[c]);
977 T_(c) = householder_left (c);
980 for (r = c + 1; r <
N; r++) {
981 if ((s = nPvt[r]) > 0) {
983 nr_double_t
t =
norm (
A_(c, r) / s);
985 y = s *
sqrt (1 - t);
986 if (fabs (y / s) < NR_TINY)
987 nPvt[
r] = euclidian_c (r, c + 1);
997 template <
class nr_type_t>
1003 for (c = 0; c <
N - 1; c++) {
1005 for (f = 0, r = c; r <
N; r++) f +=
cond_conj (
A_(r, c)) *
B_(r);
1007 for (r = c; r <
N; r++)
B_(r) -= 2.0 * f *
A_(r, c);
1011 for (r = N - 1; r >= 0; 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);
1023 template <
class nr_type_t>
1029 for (c = 0; c <
N; c++) {
1035 for (r = c + 1; r <
N; r++)
B_(r) -= f *
A_(r, c);
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);
1053 template <
class nr_type_t>
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);
1068 for (c = N - 1; c >= 0; c--) {
1073 f *=
T_(c);
B_(c) -= f;
1074 for (r = c + 1; r <
N; r++)
B_(r) -= f *
A_(r, c);
1079 for (r = 0; r <
N; r++)
X_(cMap[r]) =
B_(r);
1082 #define sign_(a) (real (a) < 0 ? -1 : 1)
1090 template <
class nr_type_t>
1094 s = euclidian_c (c, c + 1);
1095 if (s == 0 &&
imag (
A_(c, c)) == 0) {
1106 for (
int r = c + 1;
r <
N;
r++)
A_(
r, c) /= b;
1116 template <
class nr_type_t>
1119 nr_type_t
t = householder_create_left (c);
1122 householder_apply_left (c, t);
1131 template <
class nr_type_t>
1134 nr_type_t
t = householder_create_right (r);
1137 householder_apply_right (r, t);
1148 template <
class nr_type_t>
1152 s = euclidian_r (r, r + 2);
1153 if (s == 0 &&
imag (
A_(r, r + 1)) == 0) {
1164 for (
int c = r + 2;
c <
N;
c++)
A_(r,
c) /= b;
1172 template <
class nr_type_t>
1177 for (r = c + 1; r <
N; r++) {
1180 for (k = c + 1; k <
N; k++) f +=
cond_conj (
A_(k, c)) *
A_(k, r);
1183 for (k = c + 1; k <
N; k++)
A_(k, r) -= f *
A_(k, c);
1189 template <
class nr_type_t>
1194 for (c = r + 1; c <
N; c++) {
1197 for (k = r + 2; k <
N; k++) f +=
cond_conj (
A_(r, k)) *
A_(c, k);
1200 for (k = r + 2; k <
N; k++)
A_(c, k) -= f *
A_(r, k);
1213 template <
class nr_type_t>
1218 for (c = r + 1; c <
N; c++) {
1221 for (k = r + 2; k <
N; k++) f +=
cond_conj (
A_(r, k)) *
V_(c, k);
1224 for (k = r + 2; k <
N; k++)
V_(c, k) -= f *
A_(r, k);
1230 template <
class nr_type_t>
1238 template <
class nr_type_t>
1241 nr_double_t Max, Min;
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;
1251 template <
class nr_type_t>
1256 for (c = 0; c <
N; c++) {
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);
1277 template <
class nr_type_t>
1290 for (i = 0; i <
N; i++) {
1291 T_(i) = householder_left (i);
1292 if (i < N - 1)
R_(i) = householder_right (i);
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));
1301 for (l = N, i = N - 1; i >= 0; l = i--) {
1303 if ((t =
R_(i)) != 0.0) {
1304 householder_apply_right_extern (i,
cond_conj (t));
1306 else for (j = l; j <
N; j++)
1307 V_(i, j) =
V_(j, i) = 0.0;
1314 for (l = N, i = N - 1; i >= 0; l = i--) {
1315 for (j = l; j <
N; j++)
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;
1321 else for (j = l; j <
N; j++)
1332 # define MAX(x,y) (((x) > (y)) ? (x) : (y))
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);
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;
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;
1371 template <
class nr_type_t>
1374 int i,
l, j, its, k,
n, MaxIters = 30;
1375 nr_double_t an, f, g,
h,
d,
c,
s, b, a;
1378 for (an = 0, i = 0; i <
N; i++)
1379 an =
MAX (an, fabs (
S_(i)) + fabs (
E_(i)));
1383 for (k = N - 1; k >= 0; k--) {
1385 for (its = 0; its <= MaxIters; its++) {
1390 for (n = 0, l = k; l >= 1; l--) {
1393 if (fabs (
E_(l)) + an == an) { split =
false;
break; }
1394 if (fabs (
S_(n)) + an == an)
break;
1402 for (i = l; i <= k; i++) {
1405 if (fabs (f) + an == an)
break;
1409 givens_apply_u (n, i, c, s);
1419 for (j = 0; j <
N; j++)
V_(k, j) = -
V_(k, j);
1423 if (its == MaxIters) {
1437 f = (b -
d) * (b + d) + (g -
h) * (g + h);
1440 f = ((a -
d) * (a + d) + h * (b / f -
h)) / a;
1446 for (j = l; j <=
n; j++) {
1464 givens_apply_v (j, i, c, s);
1482 givens_apply_u (j, i, c, s);
matrix inverse(matrix a)
Compute inverse matrix.
void ensure_diagonal_MNA(void)
matrix real(matrix a)
Real part matrix.
matrix abs(matrix a)
Computes magnitude of each matrix element.
void householder_apply_left(int, nr_type_t)
void solve_gauss_jordan(void)
void ensure_diagonal(void)
void solve_lu_crout(void)
eqnsys()
Constructor creates an unnamed instance of the eqnsys class.
nr_double_t convergence_criteria(void)
nr_type_t householder_left(int)
void factorize_lu_crout(void)
nr_type_t householder_create_left(int)
void substitute_lu_doolittle(void)
void substitute_qr_householder_ls(void)
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.
void passEquationSys(tmatrix< nr_type_t > *, tvector< nr_type_t > *, tvector< nr_type_t > *)
void householder_apply_right(int, nr_type_t)
void factorize_qr_householder(void)
void factorize_lu_doolittle(void)
nr_complex_t sign(const nr_complex_t z)
complex sign function
#define throw_exception(e)
matrix imag(matrix a)
Imaginary part matrix.
void substitute_qrh(void)
nr_complex_t sqrt(const nr_complex_t z)
Compute principal value of square root.
nr_double_t xhypot(const nr_complex_t a, const nr_complex_t b)
Euclidean distance function for complex argument.
void givens_apply_v(int, int, nr_double_t, nr_double_t)
void substitute_svd(void)
void solve_lu_doolittle(void)
~eqnsys()
Destructor deletes the eqnsys class object.
static void euclidian_update(nr_double_t a, nr_double_t &n, nr_double_t &scale)
#define Swap(type, a, b)
Little helper macro.
void preconditioner(void)
nr_double_t euclidian_r(int, int c=1)
void chop_svd(void)
Annihilates near-zero singular values.
vector diff(vector, vector, int n=1)
void householder_apply_right_extern(int, nr_type_t)
nr_type_t cond_conj(nr_type_t t)
void substitute_qr_householder(void)
nr_double_t euclidian_c(int, int r=1)
nr_double_t norm(const nr_complex_t z)
Compute euclidian norm of complex number.
void diagonalize_svd(void)
int countPairs(int, int &, int &)
void substitute_lu_crout(void)
void setText(const char *,...)
nr_type_t householder_create_right(int)
matrix conj(matrix a)
Conjugate complex matrix.
#define VIRTUAL_RES(txt, i)
void logprint(int level, const char *format,...)
void givens_apply_u(int, int, nr_double_t, nr_double_t)
nr_type_t householder_right(int)
void solve_iterative(void)
#define MAX(x, y)
Maximum of x and y.