/* nag_dtgevc (f08ykc) Example Program.
*
* Copyright 2014 Numerical Algorithms Group.
*
* Mark 7, 2001.
*/
#include <stdio.h>
#include <nag.h>
#include <nag_stdlib.h>
#include <nagf06.h>
#include <nagf08.h>
#include <nagf16.h>
#include <nagx04.h>
static Integer normalize_vectors(Nag_OrderType order, Integer n, double qz[],
double alphai[], const char* title);
int main(void)
{
/* Scalars */
Integer i, icols, ihi, ilo, irows, j, m, n, pda, pdb, pdq, pdz;
Integer exit_status = 0;
Nag_Boolean ileft, iright;
NagError fail;
Nag_OrderType order;
/* Arrays */
double *a = 0, *alphai = 0, *alphar = 0, *b = 0, *beta = 0;
double *lscale = 0, *q = 0, *rscale = 0, *tau = 0, *z = 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_dtgevc (f08ykc) Example Program Results\n\n");
/* ileft is Nag_TRUE if left eigenvectors are required */
/* iright is Nag_TRUE if right eigenvectors are required */
ileft = Nag_TRUE;
iright = Nag_TRUE;
/* Skip heading in data file */
scanf("%*[^\n] ");
scanf("%ld %*[^\n] ", &n);
pda = n;
pdb = n;
pdq = n;
pdz = n;
/* Allocate memory */
if (
!(a = NAG_ALLOC(n * n, double)) ||
!(b = NAG_ALLOC(n * n, double)) ||
!(q = NAG_ALLOC(n * n, double)) ||
!(z = NAG_ALLOC(n * n, double)) ||
!(alphai = NAG_ALLOC(n, double)) ||
!(alphar = NAG_ALLOC(n, double)) ||
!(beta = NAG_ALLOC(n, double)) ||
!(lscale = NAG_ALLOC(n, double)) ||
!(rscale = NAG_ALLOC(n, double)) ||
!(tau = 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", &A(i, j));
scanf("%*[^\n] ");
/* READ matrix B from data file */
for (i = 1; i <= n; ++i)
for (j = 1; j <= n; ++j)
scanf("%lf", &B(i, j));
scanf("%*[^\n] ");
/* Balance the real general matrix pair (A,B) using nag_dggbal (f08whc). */
nag_dggbal(order, Nag_DoBoth, n, a, pda, b, pdb, &ilo, &ihi, lscale,
rscale, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_dggbal (f08whc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Print matrices A and B after balancing using
* nag_gen_real_mat_print (x04cac).
*/
fflush(stdout);
nag_gen_real_mat_print(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n, n, a,
pda, "Matrix A after balancing", 0, &fail);
if (fail.code == NE_NOERROR) {
fflush(stdout);
nag_gen_real_mat_print(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n, n, b,
pdb, "Matrix B after balancing", 0, &fail);
}
if (fail.code != NE_NOERROR) {
printf("Error from nag_gen_real_mat_print (x04cac).\n%s\n", fail.message);
exit_status = 2;
goto END;
}
printf("\n");
/* Reduce B to triangular form using QR and multipling both sides by Q^T */
irows = ihi + 1 - ilo;
icols = n + 1 - ilo;
/* nag_dgeqrf (f08aec).
* QR factorization of real general rectangular matrix
*/
nag_dgeqrf(order, irows, icols, &B(ilo, ilo), pdb, tau, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_dgeqrf (f08aec).\n%s\n", fail.message);
exit_status = 3;
goto END;
}
/* Apply the Q to matrix A - nag_dormqr (f08agc)
* as determined by nag_dgeqrf (f08aec).
*/
nag_dormqr(order, Nag_LeftSide, Nag_Trans, irows, icols, irows,
&B(ilo, ilo), pdb, tau, &A(ilo, ilo), pda, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_dormqr (f08agc).\n%s\n", fail.message);
exit_status = 4;
goto END;
}
/* Initialize Q (if left eigenvectors are required) */
if (ileft) {
/* Q = I. */
nag_dge_load(order, n, n, 0.0, 1.0, q, pdq, &fail);
/* Copy B to Q using nag_dge_copy (f16qfc). */
nag_dge_copy(order, Nag_NoTrans, irows-1, irows-1, &B(ilo+1,ilo), pdb,
&Q(ilo+1,ilo), pdq, &fail);
/* nag_dorgqr (f08afc).
* Form all or part of orthogonal Q from QR factorization
* determined by nag_dgeqrf (f08aec) or nag_dgeqpf (f08bec)
*/
nag_dorgqr(order, irows, irows, irows, &Q(ilo, ilo), pdq, tau, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_dorgqr (f08afc).\n%s\n", fail.message);
exit_status = 5;
goto END;
}
}
/* Initialize Z (if right eigenvectors are required) */
if (iright) {
/* Z = I. */
nag_dge_load(order, n, n, 0.0, 1.0, z, pdz, &fail);
}
/* Compute the generalized Hessenberg form of (A,B) */
/* nag_dgghrd (f08wec).
* Orthogonal reduction of a pair of real general matrices
* to generalized upper Hessenberg form
*/
nag_dgghrd(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_dgghrd (f08wec).\n%s\n", fail.message);
exit_status = 6;
goto END;
}
/* Matrix A in generalized Hessenberg form */
/* nag_gen_real_mat_print (x04cac), see above. */
fflush(stdout);
nag_gen_real_mat_print(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n, n, a,
pda, "Matrix A in Hessenberg form", 0, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_gen_real_mat_print (x04cac).\n%s\n", fail.message);
exit_status = 7;
goto END;
}
printf("\n");
/* Matrix B in generalized Hessenberg form */
/* nag_gen_real_mat_print (x04cac), see above. */
fflush(stdout);
nag_gen_real_mat_print(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n, n, b,
pdb, "Matrix B in Hessenberg form", 0, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_gen_real_mat_print (x04cac).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* nag_dhgeqz (f08xec).
* Eigenvalues and generalized Schur factorization of real
* generalized upper Hessenberg form reduced from a pair of
* real general matrices.
*/
nag_dhgeqz(order, Nag_Schur, Nag_AccumulateQ, Nag_AccumulateZ, n, ilo, ihi,
a, pda, b, pdb, alphar, alphai, beta, q, pdq, z, pdz, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_dhgeqz (f08xec).\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] != 0.0) {
printf(" %4ld (%7.3f,%7.3f)\n", i+1, alphar[i]/beta[i],
alphai[i]/beta[i]);
}
else
printf(" %4ldEigenvalue is infinite\n", i+1);
}
printf("\n");
/* Compute left and right generalized eigenvectors
* of the balanced matrix - nag_dtgevc (f08ykc).
*/
nag_dtgevc(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_dtgevc (f08ykc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
if (iright) {
/* Compute right eigenvectors of the original matrix pair
* supplied tonag_dggbal (f08whc) using nag_dggbak (f08wjc).
*/
nag_dggbak(order, Nag_DoBoth, Nag_RightSide, n, ilo, ihi, lscale,
rscale, n, z, pdz, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_dggbak (f08wjc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Normalize and print the right eigenvectors */
exit_status = normalize_vectors(order, n, z, alphai, "Right eigenvectors");
}
printf("\n");
/* Compute left eigenvectors of the original matrix */
if (ileft) {
/* nag_dggbak (f08wjc), see above. */
nag_dggbak(order, Nag_DoBoth, Nag_LeftSide, n, ilo, ihi, lscale, rscale, n,
q, pdq, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_dggbak (f08wjc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Normalize the left eigenvectors */
exit_status = normalize_vectors(order, n, q, alphai, "Left eigenvectors");
}
END:
NAG_FREE(a);
NAG_FREE(b);
NAG_FREE(q);
NAG_FREE(z);
NAG_FREE(alphai);
NAG_FREE(alphar);
NAG_FREE(beta);
NAG_FREE(lscale);
NAG_FREE(rscale);
NAG_FREE(tau);
return exit_status;
}
static Integer normalize_vectors(Nag_OrderType order, Integer n, double qz[],
double alphai[], const char* title)
{
/* Real eigenvectors are scaled so that the maximum value of elements is 1.0;
* each complex eigenvector z[] is normalized so that the element of largest
* magnitude is scaled to be (1.0,0.0).
*/
double a, b, u, v, r, ri;
Integer colinc, rowinc, i, j, k, indqz, errors=0;
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++) {
if (alphai[j]>=0.0) {
if (alphai[j]==0.0) {
/* Find element of eigenvector with largest absolute value using
* nag_damax_val (f16jqc).
*/
nag_damax_val(n, &qz[indqz], rowinc, &k, &r, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_damax_val (f16jqc).\n%s\n", fail.message);
errors = 1;
goto END;
}
r = qz[indqz+k];
for (i=0; i<n*rowinc; i+=rowinc) {
qz[indqz+i] = qz[indqz+i]/r;
}
} else {
/* norm of j-th complex eigenvector using nag_dge_norm (f16rac),
* stored as two arrays of length n.
*/
k = 0;
r = -1.0;
for (i=0; i<n*rowinc;i+=rowinc) {
ri = abs(qz[indqz+i])+abs(qz[indqz+colinc+i]);
if (ri>r) {
k = i;
r = ri;
}
}
a = qz[indqz+k];
b = qz[indqz+colinc+k];
r = a*a + b*b;
for (i=0; i<n*rowinc; i+=rowinc) {
u = qz[indqz+i];
v = qz[indqz+colinc+i];
qz[indqz+i] = (u*a + v*b)/r;
qz[indqz+colinc+i] = (v*a - u*b)/r;
}
indqz += colinc;
}
indqz += colinc;
}
}
/* Print the normalized eigenvectors using
* nag_gen_real_mat_print (x04cac)
*/
fflush(stdout);
nag_gen_real_mat_print(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n, n,
qz, n, title, 0, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_gen_real_mat_print (x04cac).\n%s\n", fail.message);
errors = 1;
}
END:
return errors;
}