/* nag_dsygv (f08sac) Example Program.
 *
 * Copyright 2014 Numerical Algorithms Group.
 *
 * Mark 23, 2011.
 */

#include <stdio.h>
#include <nag.h>
#include <nag_stdlib.h>
#include <nagf07.h>
#include <nagf08.h>
#include <nagf16.h>
#include <nagx04.h>
#include <nagx02.h>

int main(void)
{
  
  /* Scalars */
  double        anorm, bnorm, eps, rcond, rcondb, t1, t2, t3;
  Integer       exit_status = 0, i, j, n, pda, pdb;
  
  /* Arrays */
  double        *a = 0, *b = 0, *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_dsygv (f08sac) Example Program Results\n\n");
  
  /* Skip heading in data file */
  scanf("%*[^\n]");
  scanf("%ld%*[^\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, double)) ||
      !(b      = NAG_ALLOC(n * n, double)) ||
      !(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", &A(i, j));
      scanf("%*[^\n]");
      for (i = 1; i <= n; ++i)
        for (j = i; j <= n; ++j) scanf("%lf", &B(i, j));
    }
  else
    {
      for (i = 1; i <= n; ++i)
        for (j = 1; j <= i; ++j) scanf("%lf", &A(i, j));
      scanf("%*[^\n] ");
      for (i = 1; i <= n; ++i)
        for (j = 1; j <= i; ++j) scanf("%lf", &B(i, j));
    }
  scanf("%*[^\n] ");
  
  /* Compute the one-norms of the symmetric matrices A and B */
  nag_dsy_norm(order, Nag_OneNorm, uplo, n, a, pda, &anorm, &fail);
  if (fail.code != NE_NOERROR)
    {
      printf("Error from nag_dsy_norm (f16rcc).\n%s\n", fail.message);
      exit_status = 1;
      goto END;
    }
  
  nag_dsy_norm(order, Nag_OneNorm, uplo, n, b, pdb, &bnorm, &fail);
  if (fail.code != NE_NOERROR)
    {
      printf("Error from nag_dsy_norm (f16rcc).\n%s\n", fail.message);
      exit_status = 1;
      goto END;
    }
  
  /* Solve the generalized symmetric eigenvalue problem  A*x = lambda*B*x 
   * using nag_dsygv (f08sac). 
   */
  nag_dsygv(order, 1, Nag_DoBoth, uplo, n, a, pda, b, pdb, w, &fail);
  if (fail.code != NE_NOERROR)
    {
      printf("Error from nag_dsygv (f08sac).\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) = A(i, j) / A(1,j);

  /* Print eigensolution */
  printf(" Eigenvalues\n   ");
  for (j = 0; j < n; ++j) printf(" %10.4f%s", w[j], j%6 == 5?"\n":"");
  printf("\n\n");

  /* 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, a,
                         pda, "Eigenvectors", 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;
    }
      
  /* Estimate the reciprocal condition number of the Cholesky factor of B.  
   * nag_dtrcon (f07tgc)
   * Note that: cond(B) = 1.0/(rcond*rcond)
   */
  nag_dtrcon(order, Nag_OneNorm, uplo, Nag_NonUnitDiag, n, b, pdb, &rcond,
             &fail);
  if (fail.code != NE_NOERROR)
    {
      printf("Error from nag_dtrcon (f07tgc).\n%s\n", fail.message);
      exit_status = 1;
      goto END;
    }
  rcondb = rcond * rcond;
  printf("\nEstimate of reciprocal condition number for B\n%15.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-lambda*B
   * 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 = eps / rcondb;
  t2 = anorm / bnorm;
  t3 = t2 / rcond;
  for (i = 0; i < n; ++i)
    {
      eerbnd[i] = t1 * (t2 + abs(w[i]));
      zerbnd[i] = t1 * (t3 + abs(w[i]))/rcondz[i];
    }

  /* Print the approximate error bounds for the eigenvalues and vectors */
  printf("\nError estimates for the eigenvalues\n    ");
  for (i = 0; i < n; ++i) printf(" %10.1e%s", eerbnd[i], i%6 == 5?"\n":"");

  printf("\n\nError estimates for the eigenvectors\n    ");
  for (i = 0; i < n; ++i) printf(" %10.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;
}