/* nag_ztgevc (f08yxc) Example Program.
*
* Copyright 2014 Numerical Algorithms Group.
*
* Mark 7, 2001.
*/
#include <stdio.h>
#include <nag.h>
#include <nag_stdlib.h>
#include <naga02.h>
#include <nagf06.h>
#include <nagf08.h>
#include <nagf16.h>
#include <nagx04.h>
#include <naga02.h>
static Integer normalize_vectors(Nag_OrderType order, Integer n, Complex qz[],
const char* title);
int main(void)
{
/* Scalars */
Integer i, icols, ihi, ilo, irows, j, m, n, pda, pdb, pdq, pdz;
Integer exit_status = 0;
Complex e, one, zero;
Nag_Boolean ileft, iright;
NagError fail;
Nag_OrderType order;
/* Arrays */
Complex *a = 0, *alpha = 0, *b = 0, *beta = 0, *q = 0, *tau = 0;
Complex *z = 0;
double *lscale = 0, *rscale = 0;
#ifdef NAG_COLUMN_MAJOR
#define A(I, J) a[(J-1)*pda + I - 1]
#define B(I, J) b[(J-1)*pdb + I - 1]
#define Q(I, J) q[(J-1)*pdq + I - 1]
#define Z(I, J) z[(J-1)*pdz + I - 1]
order = Nag_ColMajor;
#else
#define A(I, J) a[(I-1)*pda + J - 1]
#define B(I, J) b[(I-1)*pdb + J - 1]
#define Q(I, J) q[(I-1)*pdq + J - 1]
#define Z(I, J) z[(I-1)*pdz + J - 1]
order = Nag_RowMajor;
#endif
INIT_FAIL(fail);
printf("nag_ztgevc (f08yxc) Example Program Results\n\n");
/* ileft is true if left eigenvectors are required;
* iright is true if right eigenvectors are required.
*/
ileft = Nag_TRUE;
iright = Nag_TRUE;
zero = nag_complex(0.0,0.0);
one = nag_complex(1.0,0.0);
/* Skip heading in data file and read matrix size.*/
scanf("%*[^\n] ");
scanf("%ld%*[^\n] ", &n);
pda = n;
pdb = n;
pdq = n;
pdz = n;
/* Allocate memory */
if (
!(a = NAG_ALLOC(n * n, Complex)) ||
!(b = NAG_ALLOC(n * n, Complex)) ||
!(q = NAG_ALLOC(n * n, Complex)) ||
!(z = NAG_ALLOC(n * n, Complex)) ||
!(alpha = NAG_ALLOC(n, Complex)) ||
!(beta = NAG_ALLOC(n, Complex)) ||
!(tau = NAG_ALLOC(n, Complex)) ||
!(lscale = NAG_ALLOC(n, double)) ||
!(rscale = NAG_ALLOC(n, double)))
{
printf("Allocation failure\n");
exit_status = -1;
goto END;
}
/* READ matrix A from data file */
for (i = 1; i <= n; ++i)
for (j = 1; j <= n; ++j)
scanf(" ( %lf , %lf )", &A(i, j).re, &A(i, j).im);
scanf("%*[^\n] ");
/* READ matrix B from data file */
for (i = 1; i <= n; ++i)
for (j = 1; j <= n; ++j)
scanf(" ( %lf , %lf )", &B(i, j).re, &B(i, j).im);
scanf("%*[^\n] ");
/* Balance pair (A,B) of complex general matrices using
* nag_zggbal (f08wvc).
*/
nag_zggbal(order, Nag_DoBoth, n, a, pda, b, pdb, &ilo, &ihi, lscale,
rscale, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zggbal (f08wvc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Print complex general matrices A and B after balancing using
* nag_gen_complx_mat_print_comp (x04dbc).
*/
fflush(stdout);
nag_gen_complx_mat_print_comp(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n,
n, a, pda, Nag_BracketForm, "%7.4f",
"Matrix A after balancing",
Nag_IntegerLabels, 0, Nag_IntegerLabels, 0, 80,
0, 0, &fail);
if (fail.code == NE_NOERROR) {
fflush(stdout);
nag_gen_complx_mat_print_comp(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n,
n, b, pdb, Nag_BracketForm, "%7.4f",
"Matrix B after balancing",
Nag_IntegerLabels, 0, Nag_IntegerLabels, 0,
80, 0, 0, &fail);
}
if (fail.code != NE_NOERROR) {
printf("Error from nag_gen_complx_mat_print_comp (x04dbc).\n%s\n",
fail.message);
exit_status = 1;
goto END;
}
printf("\n");
/* Reduce B to triangular form using QR and premultiply A by Q^H. */
irows = ihi + 1 - ilo;
icols = n + 1 - ilo;
/* nag_zgeqrf (f08asc).
* QR factorization of complex general rectangular matrix B.
*/
nag_zgeqrf(order, irows, icols, &B(ilo, ilo), pdb, tau, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zgeqrf (f08asc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Apply the orthogonal transformation Q^H to matrix A using
* nag_zunmqr (f08auc).
*/
nag_zunmqr(order, Nag_LeftSide, Nag_ConjTrans, irows, icols, irows,
&B(ilo, ilo), pdb, tau, &A(ilo, ilo), pda, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zunmqr (f08auc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Initialize Q (if left eigenvectors are required) */
if (ileft) {
/* Q = I */
nag_zge_load(order, n, n, zero, one, q, pdq, &fail);
/* Q = B using nag_zge_copy (f16tfc). */
nag_zge_copy(order, Nag_NoTrans, irows-1, irows-1, &B(ilo+1,ilo), pdb,
&Q(ilo+1,ilo), pdq, &fail);
/* Form Q from QR factorization using nag_zungqr (f08atc). */
nag_zungqr(order, irows, irows, irows, &Q(ilo, ilo), pdq, tau, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zungqr (f08atc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
}
if (iright) {
/* Z = I. */
nag_zge_load(order, n, n, zero, one, z, pdz, &fail);
}
/* Compute the generalized Hessenberg form of (A,B) by Unitary reduction
* using nag_zgghrd (f08wsc).
*/
nag_zgghrd(order, Nag_UpdateSchur, Nag_UpdateZ, n, ilo, ihi, a, pda, b, pdb,
q, pdq, z, pdz, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zgghrd (f08wsc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Print generalized Hessenberg form of (A,B) using
* nag_gen_complx_mat_print_comp (x04dbc).
*/
fflush(stdout);
nag_gen_complx_mat_print_comp(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n,
n, a, pda, Nag_BracketForm, "%7.3f",
"Matrix A in Hessenberg form",
Nag_IntegerLabels, 0, Nag_IntegerLabels, 0, 80,
0, 0, &fail);
if (fail.code == NE_NOERROR) {
printf("\n");
fflush(stdout);
nag_gen_complx_mat_print_comp(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n,
n, b, pdb, Nag_BracketForm, "%7.3f",
"Matrix B in Hessenberg form",
Nag_IntegerLabels, 0, Nag_IntegerLabels, 0,
80, 0, 0, &fail);
}
if (fail.code != NE_NOERROR) {
printf("Error from nag_gen_complx_mat_print_comp (x04dbc).\n%s\n",
fail.message);
exit_status = 1;
goto END;
}
/* Compute the generalized Schur form - nag_zhgeqz (f08xsc).
* Eigenvalues and generalized Schur factorization of
* complex generalized upper Hessenberg form reduced from a
* pair of complex general matrices
*/
nag_zhgeqz(order, Nag_Schur, Nag_AccumulateQ, Nag_AccumulateZ, n, ilo, ihi,
a, pda, b, pdb, alpha, beta, q, pdq, z, pdz, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zhgeqz (f08xsc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Print the generalized eigenvalue parameters */
printf("\n Generalized eigenvalues\n");
for (i = 0; i < n; ++i) {
if (beta[i].re != 0.0 || beta[i].im != 0.0) {
/* nag_complex_divide (a02cdc) - Quotient of two complex numbers. */
e = nag_complex_divide(alpha[i], beta[i]);
printf(" %4ld (%7.3f,%7.3f)\n", i+1, e.re, e.im);
}
else
printf(" %4ldEigenvalue is infinite\n", i+1);
}
printf("\n");
/* nag_ztgevc (f08yxc).
* Left and right eigenvectors of a pair of complex upper
* triangular matrices
*/
nag_ztgevc(order, Nag_BothSides, Nag_BackTransform, NULL, n, a, pda,
b, pdb, q, pdq, z, pdz, n, &m, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_ztgevc (f08yxc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
if (iright) {
/* nag_zggbak (f08wwc).
* Transform eigenvectors of a pair of complex balanced
* matrices to those of original matrix pair supplied to
* nag_zggbal (f08wvc)
*/
nag_zggbak(order, Nag_DoBoth, Nag_RightSide, n, ilo, ihi, lscale,
rscale, n, z, pdz, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zggbak (f08wwc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Normalize and print the right eigenvectors */
exit_status = normalize_vectors(order, n, z, "Right eigenvectors");
printf("\n");
}
/* Compute left eigenvectors of the original matrix */
if (ileft) {
/* nag_zggbak (f08wwc), see above. */
nag_zggbak(order, Nag_DoBoth, Nag_LeftSide, n, ilo, ihi, lscale,
rscale, n, q, pdq, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zggbak (f08wwc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Normalize and print the left eigenvectors */
exit_status = normalize_vectors(order, n, q, "Left eigenvectors");
}
END:
NAG_FREE(a);
NAG_FREE(b);
NAG_FREE(q);
NAG_FREE(z);
NAG_FREE(alpha);
NAG_FREE(beta);
NAG_FREE(tau);
NAG_FREE(lscale);
NAG_FREE(rscale);
return exit_status;
}
static Integer normalize_vectors(Nag_OrderType order, Integer n, Complex qz[],
const char* title)
{
/* Each complex eigenvector z[] is normalized so that the element of largest
* magnitude is scaled to be (1.0,0.0).
*/
double r;
Integer colinc, rowinc, j, k, indqz, errors=0;
Complex alpha, beta, x[1];
NagError fail;
INIT_FAIL(fail);
if (order==Nag_ColMajor) {
rowinc = 1;
colinc = n;
} else {
rowinc = n;
colinc = 1;
}
indqz = 0;
for (j=0; j<n; j++) {
/* Find element of eigenvector with largest absolute value using
* nag_zamax_val (f16jsc).
*/
nag_zamax_val(n, &qz[indqz], rowinc, &k, &r, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zamax_val (f16jac).\n%s\n", fail.message);
errors = 1;
goto END;
}
/* Use nag_complex_divide (a02cdc) to form reciprocal of qz[indqz+k]. */
beta = nag_complex_divide(nag_complex(1.0,0.0),qz[indqz+k]);
/* nag_zaxpby (f16gcc) performs y := alpha*x + beta*y;
* here to make largest element (1,0).
*/
alpha = nag_complex(0.0,0.0);
nag_zaxpby(n, alpha, x, 1, beta, &qz[indqz], rowinc, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zaxpby (f16gcc).\n%s\n", fail.message);
errors = 2;
goto END;
}
indqz += colinc;
}
/* Print the normalized eigenvectors using
* nag_gen_complx_mat_print_comp (x04dbc)
*/
fflush(stdout);
nag_gen_complx_mat_print_comp(order, Nag_GeneralMatrix, Nag_NonUnitDiag,
n, n, qz, n, Nag_BracketForm, "%7.4f",
title, Nag_IntegerLabels, 0,
Nag_IntegerLabels, 0, 80, 0, 0,
&fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_gen_complx_mat_print_comp (x04dbc).\n%s\n",
fail.message);
errors = 3;
}
END:
return errors;
}