/* nag_sparse_herm_basic_setup (f11grc) Example Program.
 *
 * Copyright 2017 Numerical Algorithms Group.
 *
 * Mark 26.1, 2017.
 */

#include <nag.h>
#include <nag_stdlib.h>
#include <naga02.h>
#include <nagf11.h>

int main(void)
{
  /* Scalars */
  Integer exit_status = 0;
  double anorm, dscale, dtol, sigerr, sigmax, sigtol, stplhs, stprhs, tol;
  Integer i, irevcm, iterm, itn, its, la, lfill, lwork,
         lwreq, maxitn, maxits, monit, n, nnz, nnzc, npivm;
  /* Arrays */
  char nag_enum_arg[100];
  Complex *a = 0, *b = 0, *work = 0, *x = 0;
  double *wgt = 0;
  Integer *icol = 0, *ipiv = 0, *irow = 0, *istr = 0;
  /* NAG types */
  Nag_SparseSym_Method method;
  Nag_SparseSym_PrecType precon;
  Nag_SparseSym_Bisection sigcmp;
  Nag_NormType norm;
  Nag_SparseSym_Weight weight;
  Nag_SparseSym_Fact mic;
  Nag_SparseSym_Piv pstrat;
  NagError fail, fail1;

  INIT_FAIL(fail);
  INIT_FAIL(fail1);

  printf("nag_sparse_herm_basic_setup (f11grc) Example Program Results\n");
  /* Skip heading in data file */
  scanf("%*[^\n]");
  scanf("%" NAG_IFMT "%*[^\n]", &n);
  scanf("%" NAG_IFMT "%*[^\n]", &nnz);
  la = 2 * nnz;
  lwork = 200;
  if (!(a = NAG_ALLOC(la, Complex)) ||
      !(b = NAG_ALLOC(n, Complex)) ||
      !(work = NAG_ALLOC(lwork, Complex)) ||
      !(x = NAG_ALLOC(n, Complex)) ||
      !(wgt = NAG_ALLOC(n, double)) ||
      !(icol = NAG_ALLOC(la, Integer)) ||
      !(ipiv = NAG_ALLOC(n, Integer)) ||
      !(irow = NAG_ALLOC(la, Integer)) || !(istr = NAG_ALLOC(n + 1, Integer))
         )
  {
    printf("Allocation failure\n");
    exit_status = -1;
    goto END;
  }
  /* Read or initialize the parameters for the iterative solver */
  /* nag_enum_name_to_value (x04nac).
   * Converts NAG enum member name to value
   */
  scanf("%99s%*[^\n] ", nag_enum_arg);
  method = (Nag_SparseSym_Method) nag_enum_name_to_value(nag_enum_arg);
  scanf("%99s%*[^\n] ", nag_enum_arg);
  precon = (Nag_SparseSym_PrecType) nag_enum_name_to_value(nag_enum_arg);
  scanf("%99s%*[^\n] ", nag_enum_arg);
  sigcmp = (Nag_SparseSym_Bisection) nag_enum_name_to_value(nag_enum_arg);
  scanf("%99s%*[^\n] ", nag_enum_arg);
  norm = (Nag_NormType) nag_enum_name_to_value(nag_enum_arg);
  scanf("%99s%*[^\n] ", nag_enum_arg);
  weight = (Nag_SparseSym_Weight) nag_enum_name_to_value(nag_enum_arg);
  scanf("%" NAG_IFMT "%*[^\n] ", &iterm);
  scanf("%lf%" NAG_IFMT "%*[^\n] ", &tol, &maxitn);
  scanf("%" NAG_IFMT "%*[^\n] ", &monit);

  /* Read the parameters for the preconditioner */
  scanf("%" NAG_IFMT "%lf%*[^\n] ", &lfill, &dtol);
  scanf("%99s%*[^\n] ", nag_enum_arg);
  mic = (Nag_SparseSym_Fact) nag_enum_name_to_value(nag_enum_arg);
  scanf("%lf%*[^\n]", &dscale);
  scanf("%99s%*[^\n] ", nag_enum_arg);
  pstrat = (Nag_SparseSym_Piv) nag_enum_name_to_value(nag_enum_arg);

  /* Read the nonzero elements of the matrix A */
  for (i = 0; i < nnz; i++)
    scanf(" ( %lf , %lf ) %" NAG_IFMT "%" NAG_IFMT "%*[^\n]",
          &a[i].re, &a[i].im, &irow[i], &icol[i]);

  /* Read right-hand side vector b and initial approximate solution x */
  for (i = 0; i < n; i++)
    scanf(" ( %lf , %lf ) ", &b[i].re, &b[i].im);
  scanf("%*[^\n]");
  for (i = 0; i < n; i++)
    scanf(" ( %lf , %lf ) ", &x[i].re, &x[i].im);

  /* Calculate incomplete Cholesky factorization preconditioner using 
   * nag_sparse_herm_chol_fac (f11jnc).
   * Complex sparse Hermitian matrix, incomplete Cholesky factorization
   */
  nag_sparse_herm_chol_fac(n, nnz, a, la, irow, icol, lfill, dtol, mic,
                           dscale, pstrat, ipiv, istr, &nnzc, &npivm, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_herm_chol_fac (f11jnc)\n%s\n",
           fail.message);
    exit_status = 1;
    goto END;
  }

  /* Initialize the solver using nag_sparse_herm_basic_setup (f11grc). */
  anorm = 0.0;
  sigmax = 0.0;
  sigtol = 0.01;
  maxits = n;
  nag_sparse_herm_basic_setup(method, precon, sigcmp, norm, weight, iterm,
                              n, tol, maxitn, anorm, sigmax, sigtol, maxits,
                              monit, &lwreq, work, lwork, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_herm_basic_setup (f11grc)\n%s\n",
           fail.message);
    exit_status = 2;
    goto END;
  }

  /* Call solver repeatedly to solve the equations
   * Note that the arrays b and x are overwritten on final exit:
   * x will contain the solution and b the residual vector
   */
  irevcm = 0;
  lwreq = lwork;
  /* First call to nag_sparse_herm_basic_solver (f11gsc).
   * Complex sparse Hermitian linear systems, preconditioned conjugate
   * gradient or Lanczos
   */
  nag_sparse_herm_basic_solver(&irevcm, x, b, wgt, work, lwreq, &fail);
  while (irevcm != 4) {
    switch (irevcm) {
    case 1:
      /* nag_sparse_herm_matvec (f11xsc).
       * Complex sparse Hermitian matrix vector multiply
       */
      nag_sparse_herm_matvec(n, nnz, a, irow, icol, Nag_SparseSym_NoCheck,
                             x, b, &fail1);
      break;
    case 2:
      /* nag_sparse_herm_precon_ichol_solve (f11jpc).
       * Solution of complex linear system involving incomplete Cholesky
       * preconditioning matrix generated by f11jnc
       */
      nag_sparse_herm_precon_ichol_solve(n, a, la, irow, icol, ipiv, istr,
                                         Nag_SparseSym_NoCheck, x, b, &fail1);
      break;
    case 3:
      /* nag_sparse_herm_basic_diagnostic (f11gtc).
       * Complex sparse Hermitian linear systems, diagnostic for f11gsc
       */
      nag_sparse_herm_basic_diagnostic(&itn, &stplhs, &stprhs, &anorm,
                                       &sigmax, &its, &sigerr, work, lwreq,
                                       &fail1);
      printf("\nMonitoring at iteration no.%4" NAG_IFMT "\n", itn);
      printf("residual norm:%14.4e\n", stplhs);
      printf("%6sSolution vector%18sResidual vector\n", "", "");
      for (i = 0; i < n; i++)
        printf(" (%13.4e, %13.4e)    (%13.4e, %13.4e)\n",
               x[i].re, x[i].im, b[i].re, b[i].im);
      printf("\n");
    }
    if (fail1.code != NE_NOERROR)
      irevcm = 6;
    /* Next call to nag_sparse_herm_basic_solver (f11gsc). */
    nag_sparse_herm_basic_solver(&irevcm, x, b, wgt, work, lwreq, &fail);
  }
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_herm_basic_solver (f11gsc).\n%s\n",
           fail.message);
    exit_status = 3;
    goto END;
  }

  /* Obtain information about the computation using
     nag_sparse_herm_basic_diagnostic (f11gtc)
   */
  nag_sparse_herm_basic_diagnostic(&itn, &stplhs, &stprhs, &anorm, &sigmax,
                                   &its, &sigerr, work, lwreq, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_herm_basic_diagnostic (f11gtc)\n%s\n",
           fail.message);
    exit_status = 4;
    goto END;
  }

  /* Print the output data */
  printf("Final Results\n");
  printf("Number of iterations for convergence:     %5" NAG_IFMT " \n", itn);
  printf("Residual norm:                            %14.4e\n", stplhs);
  printf("Right-hand side of termination criterion: %14.4e \n", stprhs);
  printf("1-norm of matrix A:                       %14.4e \n", anorm);
  printf("Largest singular value of A_bar:          %14.4e \n", sigmax);
  /* Output x */
  printf("\n%6sSolution vector%18sResidual vector\n", "", "");
  for (i = 0; i < n; i++)
    printf(" (%13.4e, %13.4e)   (%13.4e, %13.4e)\n",
           x[i].re, x[i].im, b[i].re, b[i].im);
  printf("\n");
END:
  NAG_FREE(a);
  NAG_FREE(b);
  NAG_FREE(work);
  NAG_FREE(x);
  NAG_FREE(wgt);
  NAG_FREE(icol);
  NAG_FREE(ipiv);
  NAG_FREE(irow);
  NAG_FREE(istr);
  return exit_status;
}