/* nag_zhegvd (f08sqc) Example Program.
*
* NAGPRODCODE Version.
*
* Copyright 2016 Numerical Algorithms Group.
*
* Mark 26, 2016.
*/
#include <math.h>
#include <stdio.h>
#include <nag.h>
#include <nagx04.h>
#include <nag_stdlib.h>
#include <nagf07.h>
#include <nagf08.h>
#include <nagf16.h>
#include <nagx02.h>
#include <naga02.h>
int main(void)
{
/* Scalars */
double anorm, bnorm, eps, rcond, rcondb, t1, t2, t3;
Integer i, j, n, pda, pdb;
Integer exit_status = 0;
/* Arrays */
Complex *a = 0, *b = 0;
double *eerbnd = 0, *rcondz = 0, *w = 0, *zerbnd = 0;
char nag_enum_arg[40];
/* Nag Types */
NagError fail;
Nag_OrderType order;
Nag_UploType uplo;
#ifdef NAG_COLUMN_MAJOR
#define A(I, J) a[(J-1)*pda + I - 1]
#define B(I, J) b[(J-1)*pdb + 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]
order = Nag_RowMajor;
#endif
INIT_FAIL(fail);
printf("nag_zhegvd (f08sqc) Example Program Results\n\n");
/* Skip heading in data file */
scanf("%*[^\n]");
scanf("%" NAG_IFMT "%*[^\n]", &n);
if (n < 0) {
printf("Invalid n\n");
exit_status = 1;
goto END;;
}
scanf(" %39s%*[^\n]", nag_enum_arg);
/* nag_enum_name_to_value (x04nac).
* Converts NAG enum member name to value
*/
uplo = (Nag_UploType) nag_enum_name_to_value(nag_enum_arg);
pda = n;
pdb = n;
/* Allocate memory */
if (!(a = NAG_ALLOC(n * n, Complex)) ||
!(b = NAG_ALLOC(n * n, Complex)) ||
!(eerbnd = NAG_ALLOC(n, double)) ||
!(rcondz = NAG_ALLOC(n, double)) ||
!(w = NAG_ALLOC(n, double)) || !(zerbnd = NAG_ALLOC(n, double)))
{
printf("Allocation failure\n");
exit_status = -1;
goto END;
}
/* Read the triangular parts of the matrices A and B */
if (uplo == Nag_Upper) {
for (i = 1; i <= n; ++i)
for (j = i; j <= n; ++j)
scanf(" ( %lf , %lf ) ", &A(i, j).re, &A(i, j).im);
scanf("%*[^\n]");
for (i = 1; i <= n; ++i)
for (j = i; j <= n; ++j)
scanf(" ( %lf , %lf ) ", &B(i, j).re, &B(i, j).im);
}
else {
for (i = 1; i <= n; ++i)
for (j = 1; j <= i; ++j)
scanf(" ( %lf , %lf ) ", &A(i, j).re, &A(i, j).im);
scanf("%*[^\n]");
for (i = 1; i <= n; ++i)
for (j = 1; j <= i; ++j)
scanf(" ( %lf , %lf ) ", &B(i, j).re, &B(i, j).im);
}
scanf("%*[^\n]");
/* Compute the one-norms of the symmetric matrices A and B
* using nag_zhe_norm (f16ucc).
*/
nag_zhe_norm(order, Nag_OneNorm, uplo, n, a, pda, &anorm, &fail);
nag_zhe_norm(order, Nag_OneNorm, uplo, n, b, pdb, &bnorm, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zhe_norm (f16ucc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Solve the generalized Hermitian eigenvalue problem A*B*x = lambda*x
* using nag_zhegvd (f08sqc).
*/
nag_zhegvd(order, 2, Nag_DoBoth, uplo, n, a, pda, b, pdb, w, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_zhegvd (f08sqc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Normalize the eigenvectors */
for (j = 1; j <= n; j++)
for (i = n; i >= 1; i--)
A(i, j) = nag_complex_divide(A(i, j), A(1, j));
/* Print eigensolution */
printf(" Eigenvalues\n ");
for (j = 0; j < n; ++j)
printf(" %11.4f%s", w[j], j % 6 == 5 ? "\n" : "");
printf("\n");
/* Prnit normalized vectors using nag_gen_complx_mat_print (x04dac). */
fflush(stdout);
nag_gen_complx_mat_print(order, Nag_GeneralMatrix, Nag_NonUnitDiag, n, n,
a, pda, "Eigenvectors", 0, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_gen_complx_mat_print (x04dac).\n%s\n",
fail.message);
exit_status = 1;
goto END;
}
/* Estimate the reciprocal condition number of the Cholesky factor of B.
* to estimate the reciprocal condition.
* nag_ztrcon (f07tuc)
* Note that: cond(B) = 1/rcond**2
*/
nag_ztrcon(order, Nag_OneNorm, uplo, Nag_NonUnitDiag, n, b, pdb,
&rcond, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_ztrcon (f07tuc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Print the reciprocal condition number of B */
rcondb = rcond * rcond;
printf("\nEstimate of reciprocal condition number for B\n %11.1e\n",
rcondb);
/* Get the machine precision, using nag_machine_precision (x02ajc) */
eps = nag_machine_precision;
if (rcond < eps) {
printf("\nB is very ill-conditioned, error estimates have not been "
"computed\n");
goto END;
}
/* Estimate reciprocal condition numbers for the eigenvectors of A*B-lambda*I
* nag_ddisna (f08flc)
*/
nag_ddisna(Nag_EigVecs, n, n, w, rcondz, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_ddisna (f08flc).\n%s\n", fail.message);
exit_status = 1;
goto END;
}
/* Compute the error estimates for the eigenvalues and eigenvectors. */
t1 = 1.0 / rcond;
t2 = eps * t1;
t3 = anorm * bnorm;
for (i = 0; i < n; ++i) {
eerbnd[i] = eps * (t3 + fabs(w[i]) / rcondb);
zerbnd[i] = t2 * (t3 / rcondz[i] + t1);
}
/* Print the approximate error bounds for the eigenvalues and vectors. */
printf("\nError estimates for the eigenvalues\n ");
for (i = 0; i < n; ++i)
printf(" %11.1e%s", eerbnd[i], i % 6 == 5 ? "\n" : "");
printf("\n\nError estimates for the eigenvectors\n ");
for (i = 0; i < n; ++i)
printf(" %11.1e%s", zerbnd[i], i % 6 == 5 ? "\n" : "");
printf("\n");
END:
NAG_FREE(a);
NAG_FREE(b);
NAG_FREE(eerbnd);
NAG_FREE(rcondz);
NAG_FREE(w);
NAG_FREE(zerbnd);
return exit_status;
}