/* nag_real_symm_sparse_eigensystem_iter (f12fbc) Example Program.
 *
 * Copyright 2014 Numerical Algorithms Group.
 *
 * Mark 8, 2005.
 */

#include <math.h>
#include <nag.h>
#include <nag_stdlib.h>
#include <stdio.h>
#include <nagf12.h>
#include <nagf16.h>

static void my_dgttrf(Integer, double *, double *, double *,
                      double *, Integer *, Integer *);
static void my_dgttrs(Integer, double *, double *, double *,
                      double *, Integer *, double *, double *);
static void av(Integer, Integer, double *, double *);
static void atv(Integer, Integer, double *, double *);

static int ex1(void), ex2(void);

int main(void)
{
  Integer  exit_status_ex1 = 0;
  Integer  exit_status_ex2 = 0;

  printf("nag_real_symm_sparse_eigensystem_iter (f12fbc) Example "
          "Program Results\n");
  exit_status_ex1 = ex1();
  exit_status_ex2 = ex2();

  return (exit_status_ex1 == 0 && exit_status_ex2 == 0) ? 0 : 1;
}

int ex1(void)
{
  /* Constants */
  Integer  licomm = 140, imon = 0;
  /* Scalars */
  double   estnrm, h2, sigma;
  Integer  exit_status = 0, info, irevcm, j, lcomm, n, nconv, ncv;
  Integer  nev, niter, nshift;
  /* Nag types */
  NagError fail;
  /* Arrays */
  double   *dd = 0, *dl = 0, *du = 0, *du2 = 0, *comm = 0, *eigest = 0;
  double   *eigv = 0, *resid = 0, *v = 0;
  Integer  *icomm = 0, *ipiv = 0;
  /* Pointers */
  double   *mx = 0, *x = 0, *y = 0;

  INIT_FAIL(fail);

  printf("\nExample 1\n");
  /* Skip heading in data file */
  scanf("%*[^\n] ");
  scanf("%*[^\n] ");

  /* Read values for nx, nev and cnv from data file. */
  scanf("%ld%ld%ld%*[^\n] ", &n, &nev, &ncv);

  /* Allocate memory */
  lcomm = 3*n + ncv*ncv + 8*ncv + 60;
  if (!(dd = NAG_ALLOC(n, double)) ||
      !(dl = NAG_ALLOC(n, double)) ||
      !(du = NAG_ALLOC(n, double)) ||
      !(du2 = NAG_ALLOC(n, double)) ||
      !(comm = NAG_ALLOC(lcomm, double)) ||
      !(eigv = NAG_ALLOC(ncv, double)) ||
      !(eigest = NAG_ALLOC(ncv, double)) ||
      !(resid = NAG_ALLOC(n, double)) ||
      !(v = NAG_ALLOC(n * ncv, double)) ||
      !(icomm = NAG_ALLOC(licomm, Integer)) ||
      !(ipiv = NAG_ALLOC(n, Integer)))
    {
      printf("Allocation failure\n");
      exit_status = -1;
      goto END;
    }
  /* Initialise communication arrays for problem using
     nag_real_symm_sparse_eigensystem_init (f12fac). */
  nag_real_symm_sparse_eigensystem_init(n, nev, ncv, icomm, licomm, comm,
                                        lcomm, &fail);
  if (fail.code != NE_NOERROR)
    {
      printf("Error from nag_real_symm_sparse_eigensystem_init "
              "(f12fac).\n%s\n", fail.message);
      exit_status = 1;
      goto END;
    }
  /* Select the required spectrum using
     nag_real_symm_sparse_eigensystem_option (f12fdc). */
  nag_real_symm_sparse_eigensystem_option("largest magnitude", icomm, comm,
                                          &fail);
  if (fail.code != NE_NOERROR)
    {
      printf("Error from nag_real_symm_sparse_eigensystem_option "
              "(f12fdc).\n%s\n", fail.message);
      exit_status = 1;
      goto END;
    }
  /* Select the required mode */
  nag_real_symm_sparse_eigensystem_option("shifted inverse", icomm, comm,
                                          &fail);
  if (fail.code != NE_NOERROR)
    {
      printf("Error from nag_real_symm_sparse_eigensystem_option "
              "(f12fdc).\n%s\n", fail.message);
      exit_status = 1;
      goto END;
    }

  h2 = 1.0 / (double)((n + 1) * (n + 1));
  sigma = 0.0;
  for (j = 0; j <= n-1; ++j)
    {
      dd[j] = 2.0 / h2 - sigma;
      dl[j] = -1.0 / h2;
      du[j] = dl[j];
    }
  my_dgttrf(n, dl, dd, du, du2, ipiv, &info);

  irevcm = 0;
 REVCOMLOOP:
  /* Repeated calls to reverse communication routine
     nag_real_symm_sparse_eigensystem_iter (f12fbc). */
  nag_real_symm_sparse_eigensystem_iter(&irevcm, resid, v, &x, &y, &mx,
                                        &nshift, comm, icomm, &fail);
  if (irevcm != 5)
    {
      if (irevcm == -1 || irevcm == 1)
        {
          /* Perform  y <--- OP*x = inv[A-SIGMA*I]*x. */
          /* Use my_dgttrs, a cut down C version of Lapack's dgttrs. */
          my_dgttrs(n, dl, dd, du, du2, ipiv, x, y);
        }
      else if (irevcm == 4 && imon == 1)
        {
          /* If imon=1, get monitoring information using
             nag_real_symm_sparse_eigensystem_monit (f12fec). */
          nag_real_symm_sparse_eigensystem_monit(&niter, &nconv, eigv, eigest,
                                                 icomm, comm);
          /* Compute 2-norm of Ritz estimates using
             nag_dge_norm (f16rac).*/
          nag_dge_norm(Nag_ColMajor, Nag_FrobeniusNorm, nev, 1, eigest, nev,
                       &estnrm, &fail);
          printf("Iteration %3ld, ", niter);
          printf(" No. converged = %3ld,", nconv);
          printf(" norm of estimates = %17.8e\n", estnrm);
        }
      goto REVCOMLOOP;
    }
  if (fail.code == NE_NOERROR)
    {
      /* Post-Process using nag_real_symm_sparse_eigensystem_sol
         (f12fcc) to compute eigenvalues/vectors. */
      nag_real_symm_sparse_eigensystem_sol(&nconv, eigv, v, sigma, resid, v,
                                           comm, icomm, &fail);
      printf("\n The %4ld Ritz values", nconv);
      printf(" closest to %8.4f are:\n\n", sigma);
      for (j = 0; j <= nconv-1; ++j)
        {
          printf("%8ld%5s%12.4f\n", j+1, "", eigv[j]);
        }
    }
  else
    {
      printf(" Error from "
              "nag_real_symm_sparse_eigensystem_iter (f12fbc).\n%s\n",
              fail.message);
      exit_status = 1;
      goto END;
    }
 END:
  NAG_FREE(dd);
  NAG_FREE(dl);
  NAG_FREE(du);
  NAG_FREE(du2);
  NAG_FREE(comm);
  NAG_FREE(eigv);
  NAG_FREE(eigest);
  NAG_FREE(resid);
  NAG_FREE(v);
  NAG_FREE(icomm);
  NAG_FREE(ipiv);

  return exit_status;
}

static void my_dgttrf(Integer n, double dl[], double d[],
                      double du[], double du2[], Integer ipiv[],
                      Integer *info)
{
  /* A simple C version of the Lapack routine dgttrf with argument
     checking removed */
  /* Scalars */
  double  temp, fact;
  Integer i;
  /* Function Body */
  *info = 0;
  for (i = 0; i < n; ++i)
    {
      ipiv[i] = i;
    }
  for (i = 0; i < n - 2; ++i)
    {
      du2[i] = 0.0;
    }
  for (i = 0; i < n - 2; i++)
    {
      if (fabs(d[i]) >= fabs(dl[i]))
        {
          /* No row interchange required, eliminate dl[i]. */
          if (d[i] != 0.0)
            {
              fact = dl[i] / d[i];
              dl[i] = fact;
              d[i+1] = d[i+1] - fact * du[i];
            }
        }
      else
        {
          /* Interchange rows I and I+1, eliminate dl[I] */
          fact = d[i] / dl[i];
          d[i] = dl[i];
          dl[i] = fact;
          temp = du[i];
          du[i] = d[i+1];
          d[i+1] = temp - fact*d[i+1];
          du2[i] = du[i+1];
          du[i+1] = -fact * du[i+1];
          ipiv[i] = i + 1;
        }
    }
  if (n > 1)
    {
      i = n - 2;
      if (fabs(d[i]) >= fabs(dl[i]))
        {
          if (d[i] != 0.0)
            {
              fact = dl[i] / d[i];
              dl[i] = fact;
              d[i+1] = d[i+1] - fact * du[i];
            }
        }
      else
        {
          fact = d[i] / dl[i];
          d[i] = dl[i];
          dl[i] = fact;
          temp = du[i];
          du[i] = d[i+1];
          d[i+1] = temp - fact * d[i+1];
          ipiv[i] = i + 1;
        }
    }
  /* Check for a zero on the diagonal of U. */
  for (i = 0; i < n; ++i)
    {
      if (d[i] == 0.0)
        {
          *info = i;
          goto END;
        }
    }
 END:
  return;
}

static void my_dgttrs(Integer n, double dl[], double d[],
                      double du[], double du2[], Integer ipiv[],
                      double b[], double y[])
{
  /* A simple C version of the Lapack routine dgttrs with argument
     checking removed, the number of right-hand-sides=1, Trans='N' */
  /* Scalars */
  Integer i, ip;
  double  temp;
  /* Solve L*x = b. */
  for (i = 0; i <= n - 1; ++i)
    {
      y[i] = b[i];
    }
  for (i = 0; i < n - 1; ++i)
    {
      ip = ipiv[i];
      temp = y[i+1-ip+i] - dl[i]*y[ip];
      y[i] = y[ip];
      y[i+1] = temp;
    }
  /* Solve U*x = b. */
  y[n-1] = y[n-1] / d[n-1];
  if (n > 1)
    {
      y[n-2] = (y[n-2] - du[n-2]*y[n-1])/d[n-2];
    }
  for (i = n - 3; i >= 0; --i)
    {
      y[i] = (y[i]-du[i]*y[i+1]-du2[i]*y[i+2])/d[i];
    }
  return;
}

int ex2(void)
{
  /* Constants */
  Integer  licomm = 140;
  /* Scalars */
  double   sigma = 0, axnorm;
  Integer  exit_status = 0, irevcm, j, lcomm, m, n, nconv, ncv;
  Integer  nev, nshift;
  NagError fail;
  /* Arrays */
  double   *comm = 0, *eigv = 0, *eigest = 0;
  double   *resid = 0, *v = 0, *ax = 0;
  Integer  *icomm = 0;
  /* Ponters */
  double   *mx = 0, *x = 0, *y = 0;

  INIT_FAIL(fail);

  printf("\nExample 2\n");
  /* Skip heading in data file. */
  scanf("%*[^\n] ");

  /* Read values for m, n, nev and cnv from data file. */
  scanf("%ld%ld%ld%ld*[^\n] ",
         &m, &n, &nev, &ncv);

  /* Allocate memory */
  lcomm = 3*n + ncv*ncv + 8*ncv + 60;
  if (!(comm = NAG_ALLOC(lcomm, double)) ||
      !(eigv = NAG_ALLOC(ncv, double)) ||
      !(eigest = NAG_ALLOC(ncv, double)) ||
      !(resid = NAG_ALLOC(n, double)) ||
      !(ax = NAG_ALLOC(m, double)) ||
      !(v = NAG_ALLOC(n * ncv, double)) ||
      !(icomm = NAG_ALLOC(licomm, Integer)))
    {
      printf("Allocation failure\n");
      exit_status = -1;
      goto END;
    }
  /* Initialise communication arrays for problem using
     nag_real_symm_sparse_eigensystem_init (f12fac). */
  nag_real_symm_sparse_eigensystem_init(n, nev, ncv, icomm, licomm, comm,
                                        lcomm, &fail);
  if (fail.code != NE_NOERROR)
    {
      printf("Error from nag_real_symm_sparse_eigensystem_init "
              "(f12fac).\n%s\n", fail.message);
      exit_status = 1;
      goto END;
    }
  irevcm = 0;
 REVCOMLOOP:
  /* Repeated calls to reverse communication routine
     nag_real_symm_sparse_eigensystem_iter (f12fbc). */
  nag_real_symm_sparse_eigensystem_iter(&irevcm, resid, v, &x, &y, &mx,
                                        &nshift, comm, icomm, &fail);
  if (irevcm != 5)
    {
      if (irevcm == -1 || irevcm == 1)
        {
          /* Perform matrix vector multiplication y <--- Op*x */
          av(m, n, x, ax);
          atv(m, n, ax, y);
        }
      goto REVCOMLOOP;
    }
  if (fail.code == NE_NOERROR)
    {
      /* Post-Process using nag_real_symm_sparse_eigensystem_sol
         (f12fcc) to compute singular values/vectors. */
      nag_real_symm_sparse_eigensystem_sol(&nconv, eigv, v, sigma, resid, v,
                                           comm, icomm, &fail);
      /* Singular values (squared) are returned in eigv and the
         corresponding right singular vectors are returned in the first
         nev n-length vectors in v. */
      printf("\n The %4ld leading Singular values and", nconv);
      printf(" direct residuals are:\n\n");
      for (j = 0; j <= nconv-1; ++j)
        {
          eigv[j] = sqrt(eigv[j]);
          /* Compute the left singular vectors from the formula
             u = Av/sigma
             u should have norm 1 so divide by norm(Av). */
          av(m, n, &v[j*n], ax);
          /* Compute 2-norm of Av using nag_dge_norm (f16rac).*/
          nag_dge_norm(Nag_ColMajor, Nag_FrobeniusNorm, m, 1, ax,
                       m, &axnorm, &fail);
          resid[j] = axnorm*fabs(1.0-eigv[j]/axnorm);

          printf("%8ld%5s%12.4f%5s%12.7f\n", j+1, "", eigv[j], "",
                  resid[j]);
        }
    }
  else
    {
      printf(" Error from "
              "nag_real_symm_sparse_eigensystem_iter (f12fbc).\n%s\n",
              fail.message);
      exit_status = 1;
      goto END;
    }
 END:
  NAG_FREE(comm);
  NAG_FREE(eigv);
  NAG_FREE(eigest);
  NAG_FREE(resid);
  NAG_FREE(v);
  NAG_FREE(ax);
  NAG_FREE(icomm);

  return exit_status;
}

static void av(Integer m, Integer n, double *x, double *w)
{
  /* Computes  w <- A*x. */
  /* Local Scalars */
  double  h, k, s, t;
  Integer i, j;
  h = 1.0/(double)(m+1);
  k = 1.0/(double)(n+1);
  for (i = 0; i < m; ++i)
    {
      w[i] = 0.0;
    }
  t = 0.0;

  for (j = 0; j < n; ++j)
    {
      t = t + k;
      s = 0.0;
      for (i = 0; i < j+1; i++)
        {
          s = s + h;
          w[i] = w[i] + k*s*(t-1.0)*x[j];
        }
      for (i = j+1; i < m; ++i)
        {
          s = s + h;
          w[i] = w[i] + k*t*(s-1.0)*x[j];
        }
    }
  return;
} /* av */

static void atv(Integer m, Integer n, double *x, double *y)
{
  /* Computes  y <- A'*w. */
  /* Local Scalars */
  double  h, k, s, t;
  Integer i, j;
  h = 1.0/(double)(m+1);
  k = 1.0/(double)(n+1);
  for (i = 0; i < n; ++i)
    {
      y[i] = 0.0;
    }
  t = 0.0;
  for (j = 0; j < n; ++j)
    {
      t = t + k;
      s = 0.0;
      for (i = 0; i < j+1; ++i)
        {
          s = s + h;
          y[j] = y[j] + k*s*(t-1.0)*x[i];
        }
      for (i = j+1; i < m; ++i)
        {
          s = s + h;
          y[j] = y[j] + k*t*(s-1.0)*x[i];
        }
    }
  return;
} /* atv */