1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986
|
/* eigen/nonsymmv.c
*
* Copyright (C) 2006 Patrick Alken
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or (at
* your option) any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/
#include <config.h>
#include <stdlib.h>
#include <math.h>
#include <gsl/gsl_complex.h>
#include <gsl/gsl_complex_math.h>
#include <gsl/gsl_eigen.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_blas.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_vector_complex.h>
#include <gsl/gsl_matrix.h>
/*
* This module computes the eigenvalues and eigenvectors of a real
* nonsymmetric matrix.
*
* This file contains routines based on original code from LAPACK
* which is distributed under the modified BSD license. The LAPACK
* routines used are DTREVC and DLALN2.
*/
#define GSL_NONSYMMV_SMLNUM (2.0 * GSL_DBL_MIN)
#define GSL_NONSYMMV_BIGNUM ((1.0 - GSL_DBL_EPSILON) / GSL_NONSYMMV_SMLNUM)
static void nonsymmv_get_right_eigenvectors(gsl_matrix *T, gsl_matrix *Z,
gsl_vector_complex *eval,
gsl_matrix_complex *evec,
gsl_eigen_nonsymmv_workspace *w);
static void nonsymmv_normalize_eigenvectors(gsl_vector_complex *eval,
gsl_matrix_complex *evec);
/*
gsl_eigen_nonsymmv_alloc()
Allocate a workspace for solving the nonsymmetric eigenvalue problem.
The size of this workspace is O(5n).
Inputs: n - size of matrices
Return: pointer to workspace
*/
gsl_eigen_nonsymmv_workspace *
gsl_eigen_nonsymmv_alloc(const size_t n)
{
gsl_eigen_nonsymmv_workspace *w;
if (n == 0)
{
GSL_ERROR_NULL ("matrix dimension must be positive integer",
GSL_EINVAL);
}
w = (gsl_eigen_nonsymmv_workspace *)
calloc (1, sizeof (gsl_eigen_nonsymmv_workspace));
if (w == 0)
{
GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
}
w->size = n;
w->Z = NULL;
w->nonsymm_workspace_p = gsl_eigen_nonsymm_alloc(n);
if (w->nonsymm_workspace_p == 0)
{
gsl_eigen_nonsymmv_free(w);
GSL_ERROR_NULL ("failed to allocate space for nonsymm workspace", GSL_ENOMEM);
}
/*
* set parameters to compute the full Schur form T and balance
* the matrices
*/
gsl_eigen_nonsymm_params(1, 0, w->nonsymm_workspace_p);
w->work = gsl_vector_alloc(n);
w->work2 = gsl_vector_alloc(n);
w->work3 = gsl_vector_alloc(n);
if (w->work == 0 || w->work2 == 0 || w->work3 == 0)
{
gsl_eigen_nonsymmv_free(w);
GSL_ERROR_NULL ("failed to allocate space for nonsymmv additional workspace", GSL_ENOMEM);
}
return (w);
} /* gsl_eigen_nonsymmv_alloc() */
/*
gsl_eigen_nonsymmv_free()
Free workspace w
*/
void
gsl_eigen_nonsymmv_free (gsl_eigen_nonsymmv_workspace * w)
{
RETURN_IF_NULL (w);
if (w->nonsymm_workspace_p)
gsl_eigen_nonsymm_free(w->nonsymm_workspace_p);
if (w->work)
gsl_vector_free(w->work);
if (w->work2)
gsl_vector_free(w->work2);
if (w->work3)
gsl_vector_free(w->work3);
free(w);
} /* gsl_eigen_nonsymmv_free() */
/*
gsl_eigen_nonsymmv_params()
Set some parameters which define how we solve the eigenvalue
problem.
Inputs: balance - 1 if we want to balance the matrix, 0 if not
w - nonsymmv workspace
*/
void
gsl_eigen_nonsymmv_params (const int balance,
gsl_eigen_nonsymmv_workspace *w)
{
gsl_eigen_nonsymm_params(1, balance, w->nonsymm_workspace_p);
} /* gsl_eigen_nonsymm_params() */
/*
gsl_eigen_nonsymmv()
Solve the nonsymmetric eigensystem problem
A x = \lambda x
for the eigenvalues \lambda and right eigenvectors x
Inputs: A - general real matrix
eval - where to store eigenvalues
evec - where to store eigenvectors
w - workspace
Return: success or error
*/
int
gsl_eigen_nonsymmv (gsl_matrix * A, gsl_vector_complex * eval,
gsl_matrix_complex * evec,
gsl_eigen_nonsymmv_workspace * w)
{
const size_t N = A->size1;
/* check matrix and vector sizes */
if (N != A->size2)
{
GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
}
else if (eval->size != N)
{
GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
}
else if (evec->size1 != evec->size2)
{
GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
}
else if (evec->size1 != N)
{
GSL_ERROR ("eigenvector matrix has wrong size", GSL_EBADLEN);
}
else
{
int s;
gsl_matrix Z;
/*
* We need a place to store the Schur vectors, so we will
* treat evec as a real matrix and store them in the left
* half - the factor of 2 in the tda corresponds to the
* complex multiplicity
*/
Z.size1 = N;
Z.size2 = N;
Z.tda = 2 * N;
Z.data = evec->data;
Z.block = 0;
Z.owner = 0;
/* compute eigenvalues, Schur form, and Schur vectors */
s = gsl_eigen_nonsymm_Z(A, eval, &Z, w->nonsymm_workspace_p);
if (w->Z)
{
/*
* save the Schur vectors in user supplied matrix, since
* they will be destroyed when computing eigenvectors
*/
gsl_matrix_memcpy(w->Z, &Z);
}
/* only compute eigenvectors if we found all eigenvalues */
if (s == GSL_SUCCESS)
{
/* compute eigenvectors */
nonsymmv_get_right_eigenvectors(A, &Z, eval, evec, w);
/* normalize so that Euclidean norm is 1 */
nonsymmv_normalize_eigenvectors(eval, evec);
}
return s;
}
} /* gsl_eigen_nonsymmv() */
/*
gsl_eigen_nonsymmv_Z()
Compute eigenvalues and eigenvectors of a real nonsymmetric matrix
and also save the Schur vectors. See comments in gsl_eigen_nonsymm_Z
for more information.
Inputs: A - real nonsymmetric matrix
eval - where to store eigenvalues
evec - where to store eigenvectors
Z - where to store Schur vectors
w - nonsymmv workspace
Return: success or error
*/
int
gsl_eigen_nonsymmv_Z (gsl_matrix * A, gsl_vector_complex * eval,
gsl_matrix_complex * evec, gsl_matrix * Z,
gsl_eigen_nonsymmv_workspace * w)
{
/* check matrix and vector sizes */
if (A->size1 != A->size2)
{
GSL_ERROR ("matrix must be square to compute eigenvalues/eigenvectors", GSL_ENOTSQR);
}
else if (eval->size != A->size1)
{
GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
}
else if (evec->size1 != evec->size2)
{
GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
}
else if (evec->size1 != A->size1)
{
GSL_ERROR ("eigenvector matrix has wrong size", GSL_EBADLEN);
}
else if ((Z->size1 != Z->size2) || (Z->size1 != A->size1))
{
GSL_ERROR ("Z matrix has wrong dimensions", GSL_EBADLEN);
}
else
{
int s;
w->Z = Z;
s = gsl_eigen_nonsymmv(A, eval, evec, w);
w->Z = NULL;
return s;
}
} /* gsl_eigen_nonsymmv_Z() */
/********************************************
* INTERNAL ROUTINES *
********************************************/
/*
nonsymmv_get_right_eigenvectors()
Compute the right eigenvectors of the Schur form T and then
backtransform them using the Schur vectors to get right eigenvectors of
the original matrix.
Inputs: T - Schur form
Z - Schur vectors
eval - where to store eigenvalues (to ensure that the
correct eigenvalue is stored in the same position
as the eigenvectors)
evec - where to store eigenvectors
w - nonsymmv workspace
Return: none
Notes: 1) based on LAPACK routine DTREVC - the algorithm used is
backsubstitution on the upper quasi triangular system T
followed by backtransformation by Z to get vectors of the
original matrix.
2) The Schur vectors in Z are destroyed and replaced with
eigenvectors stored with the same storage scheme as DTREVC.
The eigenvectors are also stored in 'evec'
3) The matrix T is unchanged on output
4) Each eigenvector is normalized so that the element of
largest magnitude has magnitude 1; here the magnitude of
a complex number (x,y) is taken to be |x| + |y|
*/
static void
nonsymmv_get_right_eigenvectors(gsl_matrix *T, gsl_matrix *Z,
gsl_vector_complex *eval,
gsl_matrix_complex *evec,
gsl_eigen_nonsymmv_workspace *w)
{
const size_t N = T->size1;
const double smlnum = GSL_DBL_MIN * N / GSL_DBL_EPSILON;
const double bignum = (1.0 - GSL_DBL_EPSILON) / smlnum;
int i; /* looping */
size_t iu, /* looping */
ju,
ii;
gsl_complex lambda; /* current eigenvalue */
double lambda_re, /* Re(lambda) */
lambda_im; /* Im(lambda) */
gsl_matrix_view Tv, /* temporary views */
Zv;
gsl_vector_view y, /* temporary views */
y2,
ev,
ev2;
double dat[4], /* scratch arrays */
dat_X[4];
double scale; /* scale factor */
double xnorm; /* |X| */
gsl_vector_complex_view ecol, /* column of evec */
ecol2;
int complex_pair; /* complex eigenvalue pair? */
double smin;
/*
* Compute 1-norm of each column of upper triangular part of T
* to control overflow in triangular solver
*/
gsl_vector_set(w->work3, 0, 0.0);
for (ju = 1; ju < N; ++ju)
{
gsl_vector_set(w->work3, ju, 0.0);
for (iu = 0; iu < ju; ++iu)
{
gsl_vector_set(w->work3, ju,
gsl_vector_get(w->work3, ju) +
fabs(gsl_matrix_get(T, iu, ju)));
}
}
for (i = (int) N - 1; i >= 0; --i)
{
iu = (size_t) i;
/* get current eigenvalue and store it in lambda */
lambda_re = gsl_matrix_get(T, iu, iu);
if (iu != 0 && gsl_matrix_get(T, iu, iu - 1) != 0.0)
{
lambda_im = sqrt(fabs(gsl_matrix_get(T, iu, iu - 1))) *
sqrt(fabs(gsl_matrix_get(T, iu - 1, iu)));
}
else
{
lambda_im = 0.0;
}
GSL_SET_COMPLEX(&lambda, lambda_re, lambda_im);
smin = GSL_MAX(GSL_DBL_EPSILON * (fabs(lambda_re) + fabs(lambda_im)),
smlnum);
smin = GSL_MAX(smin, GSL_NONSYMMV_SMLNUM);
if (lambda_im == 0.0)
{
int k, l;
gsl_vector_view bv, xv;
/* real eigenvector */
/*
* The ordering of eigenvalues in 'eval' is arbitrary and
* does not necessarily follow the Schur form T, so store
* lambda in the right slot in eval to ensure it corresponds
* to the eigenvector we are about to compute
*/
gsl_vector_complex_set(eval, iu, lambda);
/*
* We need to solve the system:
*
* (T(1:iu-1, 1:iu-1) - lambda*I)*X = -T(1:iu-1,iu)
*/
/* construct right hand side */
for (k = 0; k < i; ++k)
{
gsl_vector_set(w->work,
(size_t) k,
-gsl_matrix_get(T, (size_t) k, iu));
}
gsl_vector_set(w->work, iu, 1.0);
for (l = i - 1; l >= 0; --l)
{
size_t lu = (size_t) l;
if (lu == 0)
complex_pair = 0;
else
complex_pair = gsl_matrix_get(T, lu, lu - 1) != 0.0;
if (!complex_pair)
{
double x;
/*
* 1-by-1 diagonal block - solve the system:
*
* (T_{ll} - lambda)*x = -T_{l(iu)}
*/
Tv = gsl_matrix_submatrix(T, lu, lu, 1, 1);
bv = gsl_vector_view_array(dat, 1);
gsl_vector_set(&bv.vector, 0,
gsl_vector_get(w->work, lu));
xv = gsl_vector_view_array(dat_X, 1);
gsl_schur_solve_equation(1.0,
&Tv.matrix,
lambda_re,
1.0,
1.0,
&bv.vector,
&xv.vector,
&scale,
&xnorm,
smin);
/* scale x to avoid overflow */
x = gsl_vector_get(&xv.vector, 0);
if (xnorm > 1.0)
{
if (gsl_vector_get(w->work3, lu) > bignum / xnorm)
{
x /= xnorm;
scale /= xnorm;
}
}
if (scale != 1.0)
{
gsl_vector_view wv;
wv = gsl_vector_subvector(w->work, 0, iu + 1);
gsl_blas_dscal(scale, &wv.vector);
}
gsl_vector_set(w->work, lu, x);
if (lu > 0)
{
gsl_vector_view v1, v2;
/* update right hand side */
v1 = gsl_matrix_subcolumn(T, lu, 0, lu);
v2 = gsl_vector_subvector(w->work, 0, lu);
gsl_blas_daxpy(-x, &v1.vector, &v2.vector);
} /* if (l > 0) */
} /* if (!complex_pair) */
else
{
double x11, x21;
/*
* 2-by-2 diagonal block
*/
Tv = gsl_matrix_submatrix(T, lu - 1, lu - 1, 2, 2);
bv = gsl_vector_view_array(dat, 2);
gsl_vector_set(&bv.vector, 0,
gsl_vector_get(w->work, lu - 1));
gsl_vector_set(&bv.vector, 1,
gsl_vector_get(w->work, lu));
xv = gsl_vector_view_array(dat_X, 2);
gsl_schur_solve_equation(1.0,
&Tv.matrix,
lambda_re,
1.0,
1.0,
&bv.vector,
&xv.vector,
&scale,
&xnorm,
smin);
/* scale X(1,1) and X(2,1) to avoid overflow */
x11 = gsl_vector_get(&xv.vector, 0);
x21 = gsl_vector_get(&xv.vector, 1);
if (xnorm > 1.0)
{
double beta;
beta = GSL_MAX(gsl_vector_get(w->work3, lu - 1),
gsl_vector_get(w->work3, lu));
if (beta > bignum / xnorm)
{
x11 /= xnorm;
x21 /= xnorm;
scale /= xnorm;
}
}
/* scale if necessary */
if (scale != 1.0)
{
gsl_vector_view wv;
wv = gsl_vector_subvector(w->work, 0, iu + 1);
gsl_blas_dscal(scale, &wv.vector);
}
gsl_vector_set(w->work, lu - 1, x11);
gsl_vector_set(w->work, lu, x21);
/* update right hand side */
if (lu > 1)
{
gsl_vector_view v1, v2;
v1 = gsl_matrix_subcolumn(T, lu - 1, 0, lu - 1);
v2 = gsl_vector_subvector(w->work, 0, lu - 1);
gsl_blas_daxpy(-x11, &v1.vector, &v2.vector);
v1 = gsl_matrix_subcolumn(T, lu, 0, lu - 1);
gsl_blas_daxpy(-x21, &v1.vector, &v2.vector);
}
--l;
} /* if (complex_pair) */
} /* for (l = i - 1; l >= 0; --l) */
/*
* At this point, w->work is an eigenvector of the
* Schur form T. To get an eigenvector of the original
* matrix, we multiply on the left by Z, the matrix of
* Schur vectors
*/
ecol = gsl_matrix_complex_column(evec, iu);
y = gsl_matrix_column(Z, iu);
if (iu > 0)
{
gsl_vector_view x;
Zv = gsl_matrix_submatrix(Z, 0, 0, N, iu);
x = gsl_vector_subvector(w->work, 0, iu);
/* compute Z * w->work and store it in Z(:,iu) */
gsl_blas_dgemv(CblasNoTrans,
1.0,
&Zv.matrix,
&x.vector,
gsl_vector_get(w->work, iu),
&y.vector);
} /* if (iu > 0) */
/* store eigenvector into evec */
ev = gsl_vector_complex_real(&ecol.vector);
ev2 = gsl_vector_complex_imag(&ecol.vector);
scale = 0.0;
for (ii = 0; ii < N; ++ii)
{
double a = gsl_vector_get(&y.vector, ii);
/* store real part of eigenvector */
gsl_vector_set(&ev.vector, ii, a);
/* set imaginary part to 0 */
gsl_vector_set(&ev2.vector, ii, 0.0);
if (fabs(a) > scale)
scale = fabs(a);
}
if (scale != 0.0)
scale = 1.0 / scale;
/* scale by magnitude of largest element */
gsl_blas_dscal(scale, &ev.vector);
} /* if (GSL_IMAG(lambda) == 0.0) */
else
{
gsl_vector_complex_view bv, xv;
size_t k;
int l;
gsl_complex lambda2;
/* complex eigenvector */
/*
* Store the complex conjugate eigenvalues in the right
* slots in eval
*/
GSL_SET_REAL(&lambda2, GSL_REAL(lambda));
GSL_SET_IMAG(&lambda2, -GSL_IMAG(lambda));
gsl_vector_complex_set(eval, iu - 1, lambda);
gsl_vector_complex_set(eval, iu, lambda2);
/*
* First solve:
*
* [ T(i:i+1,i:i+1) - lambda*I ] * X = 0
*/
if (fabs(gsl_matrix_get(T, iu - 1, iu)) >=
fabs(gsl_matrix_get(T, iu, iu - 1)))
{
gsl_vector_set(w->work, iu - 1, 1.0);
gsl_vector_set(w->work2, iu,
lambda_im / gsl_matrix_get(T, iu - 1, iu));
}
else
{
gsl_vector_set(w->work, iu - 1,
-lambda_im / gsl_matrix_get(T, iu, iu - 1));
gsl_vector_set(w->work2, iu, 1.0);
}
gsl_vector_set(w->work, iu, 0.0);
gsl_vector_set(w->work2, iu - 1, 0.0);
/* construct right hand side */
for (k = 0; k < iu - 1; ++k)
{
gsl_vector_set(w->work, k,
-gsl_vector_get(w->work, iu - 1) *
gsl_matrix_get(T, k, iu - 1));
gsl_vector_set(w->work2, k,
-gsl_vector_get(w->work2, iu) *
gsl_matrix_get(T, k, iu));
}
/*
* We must solve the upper quasi-triangular system:
*
* [ T(1:i-2,1:i-2) - lambda*I ] * X = s*(work + i*work2)
*/
for (l = i - 2; l >= 0; --l)
{
size_t lu = (size_t) l;
if (lu == 0)
complex_pair = 0;
else
complex_pair = gsl_matrix_get(T, lu, lu - 1) != 0.0;
if (!complex_pair)
{
gsl_complex bval;
gsl_complex x;
/*
* 1-by-1 diagonal block - solve the system:
*
* (T_{ll} - lambda)*x = work + i*work2
*/
Tv = gsl_matrix_submatrix(T, lu, lu, 1, 1);
bv = gsl_vector_complex_view_array(dat, 1);
xv = gsl_vector_complex_view_array(dat_X, 1);
GSL_SET_COMPLEX(&bval,
gsl_vector_get(w->work, lu),
gsl_vector_get(w->work2, lu));
gsl_vector_complex_set(&bv.vector, 0, bval);
gsl_schur_solve_equation_z(1.0,
&Tv.matrix,
&lambda,
1.0,
1.0,
&bv.vector,
&xv.vector,
&scale,
&xnorm,
smin);
if (xnorm > 1.0)
{
if (gsl_vector_get(w->work3, lu) > bignum / xnorm)
{
gsl_blas_zdscal(1.0/xnorm, &xv.vector);
scale /= xnorm;
}
}
/* scale if necessary */
if (scale != 1.0)
{
gsl_vector_view wv;
wv = gsl_vector_subvector(w->work, 0, iu + 1);
gsl_blas_dscal(scale, &wv.vector);
wv = gsl_vector_subvector(w->work2, 0, iu + 1);
gsl_blas_dscal(scale, &wv.vector);
}
x = gsl_vector_complex_get(&xv.vector, 0);
gsl_vector_set(w->work, lu, GSL_REAL(x));
gsl_vector_set(w->work2, lu, GSL_IMAG(x));
/* update the right hand side */
if (lu > 0)
{
gsl_vector_view v1, v2;
v1 = gsl_matrix_subcolumn(T, lu, 0, lu);
v2 = gsl_vector_subvector(w->work, 0, lu);
gsl_blas_daxpy(-GSL_REAL(x), &v1.vector, &v2.vector);
v2 = gsl_vector_subvector(w->work2, 0, lu);
gsl_blas_daxpy(-GSL_IMAG(x), &v1.vector, &v2.vector);
} /* if (lu > 0) */
} /* if (!complex_pair) */
else
{
gsl_complex b1, b2, x1, x2;
/*
* 2-by-2 diagonal block - solve the system
*/
Tv = gsl_matrix_submatrix(T, lu - 1, lu - 1, 2, 2);
bv = gsl_vector_complex_view_array(dat, 2);
xv = gsl_vector_complex_view_array(dat_X, 2);
GSL_SET_COMPLEX(&b1,
gsl_vector_get(w->work, lu - 1),
gsl_vector_get(w->work2, lu - 1));
GSL_SET_COMPLEX(&b2,
gsl_vector_get(w->work, lu),
gsl_vector_get(w->work2, lu));
gsl_vector_complex_set(&bv.vector, 0, b1);
gsl_vector_complex_set(&bv.vector, 1, b2);
gsl_schur_solve_equation_z(1.0,
&Tv.matrix,
&lambda,
1.0,
1.0,
&bv.vector,
&xv.vector,
&scale,
&xnorm,
smin);
x1 = gsl_vector_complex_get(&xv.vector, 0);
x2 = gsl_vector_complex_get(&xv.vector, 1);
if (xnorm > 1.0)
{
double beta;
beta = GSL_MAX(gsl_vector_get(w->work3, lu - 1),
gsl_vector_get(w->work3, lu));
if (beta > bignum / xnorm)
{
gsl_blas_zdscal(1.0/xnorm, &xv.vector);
scale /= xnorm;
}
}
/* scale if necessary */
if (scale != 1.0)
{
gsl_vector_view wv;
wv = gsl_vector_subvector(w->work, 0, iu + 1);
gsl_blas_dscal(scale, &wv.vector);
wv = gsl_vector_subvector(w->work2, 0, iu + 1);
gsl_blas_dscal(scale, &wv.vector);
}
gsl_vector_set(w->work, lu - 1, GSL_REAL(x1));
gsl_vector_set(w->work, lu, GSL_REAL(x2));
gsl_vector_set(w->work2, lu - 1, GSL_IMAG(x1));
gsl_vector_set(w->work2, lu, GSL_IMAG(x2));
/* update right hand side */
if (lu > 1)
{
gsl_vector_view v1, v2, v3, v4;
v1 = gsl_matrix_subcolumn(T, lu - 1, 0, lu - 1);
v4 = gsl_matrix_subcolumn(T, lu, 0, lu - 1);
v2 = gsl_vector_subvector(w->work, 0, lu - 1);
v3 = gsl_vector_subvector(w->work2, 0, lu - 1);
gsl_blas_daxpy(-GSL_REAL(x1), &v1.vector, &v2.vector);
gsl_blas_daxpy(-GSL_REAL(x2), &v4.vector, &v2.vector);
gsl_blas_daxpy(-GSL_IMAG(x1), &v1.vector, &v3.vector);
gsl_blas_daxpy(-GSL_IMAG(x2), &v4.vector, &v3.vector);
} /* if (lu > 1) */
--l;
} /* if (complex_pair) */
} /* for (l = i - 2; l >= 0; --l) */
/*
* At this point, work + i*work2 is an eigenvector
* of T - backtransform to get an eigenvector of the
* original matrix
*/
y = gsl_matrix_column(Z, iu - 1);
y2 = gsl_matrix_column(Z, iu);
if (iu > 1)
{
gsl_vector_view x;
/* compute real part of eigenvectors */
Zv = gsl_matrix_submatrix(Z, 0, 0, N, iu - 1);
x = gsl_vector_subvector(w->work, 0, iu - 1);
gsl_blas_dgemv(CblasNoTrans,
1.0,
&Zv.matrix,
&x.vector,
gsl_vector_get(w->work, iu - 1),
&y.vector);
/* now compute the imaginary part */
x = gsl_vector_subvector(w->work2, 0, iu - 1);
gsl_blas_dgemv(CblasNoTrans,
1.0,
&Zv.matrix,
&x.vector,
gsl_vector_get(w->work2, iu),
&y2.vector);
}
else
{
gsl_blas_dscal(gsl_vector_get(w->work, iu - 1), &y.vector);
gsl_blas_dscal(gsl_vector_get(w->work2, iu), &y2.vector);
}
/*
* Now store the eigenvectors into evec - the real parts
* are Z(:,iu - 1) and the imaginary parts are
* +/- Z(:,iu)
*/
/* get views of the two eigenvector slots */
ecol = gsl_matrix_complex_column(evec, iu - 1);
ecol2 = gsl_matrix_complex_column(evec, iu);
/*
* save imaginary part first as it may get overwritten
* when copying the real part due to our storage scheme
* in Z/evec
*/
ev = gsl_vector_complex_imag(&ecol.vector);
ev2 = gsl_vector_complex_imag(&ecol2.vector);
scale = 0.0;
for (ii = 0; ii < N; ++ii)
{
double a = gsl_vector_get(&y2.vector, ii);
scale = GSL_MAX(scale,
fabs(a) + fabs(gsl_vector_get(&y.vector, ii)));
gsl_vector_set(&ev.vector, ii, a);
gsl_vector_set(&ev2.vector, ii, -a);
}
/* now save the real part */
ev = gsl_vector_complex_real(&ecol.vector);
ev2 = gsl_vector_complex_real(&ecol2.vector);
for (ii = 0; ii < N; ++ii)
{
double a = gsl_vector_get(&y.vector, ii);
gsl_vector_set(&ev.vector, ii, a);
gsl_vector_set(&ev2.vector, ii, a);
}
if (scale != 0.0)
scale = 1.0 / scale;
/* scale by largest element magnitude */
gsl_blas_zdscal(scale, &ecol.vector);
gsl_blas_zdscal(scale, &ecol2.vector);
/*
* decrement i since we took care of two eigenvalues at
* the same time
*/
--i;
} /* if (GSL_IMAG(lambda) != 0.0) */
} /* for (i = (int) N - 1; i >= 0; --i) */
} /* nonsymmv_get_right_eigenvectors() */
/*
nonsymmv_normalize_eigenvectors()
Normalize eigenvectors so that their Euclidean norm is 1
Inputs: eval - eigenvalues
evec - eigenvectors
*/
static void
nonsymmv_normalize_eigenvectors(gsl_vector_complex *eval,
gsl_matrix_complex *evec)
{
const size_t N = evec->size1;
size_t i; /* looping */
gsl_complex ei;
gsl_vector_complex_view vi;
gsl_vector_view re, im;
double scale; /* scaling factor */
for (i = 0; i < N; ++i)
{
ei = gsl_vector_complex_get(eval, i);
vi = gsl_matrix_complex_column(evec, i);
re = gsl_vector_complex_real(&vi.vector);
if (GSL_IMAG(ei) == 0.0)
{
scale = 1.0 / gsl_blas_dnrm2(&re.vector);
gsl_blas_dscal(scale, &re.vector);
}
else if (GSL_IMAG(ei) > 0.0)
{
im = gsl_vector_complex_imag(&vi.vector);
scale = 1.0 / gsl_hypot(gsl_blas_dnrm2(&re.vector),
gsl_blas_dnrm2(&im.vector));
gsl_blas_zdscal(scale, &vi.vector);
vi = gsl_matrix_complex_column(evec, i + 1);
gsl_blas_zdscal(scale, &vi.vector);
}
}
} /* nonsymmv_normalize_eigenvectors() */
|