/* nag_sparse_nherm_precon_bdilu_solve (f11duc) Example Program.
 *
 * Copyright 2017 Numerical Algorithms Group.
 *
 * Mark 26.1, 2017.
 */

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

static void overlap(Integer *n, Integer *nnz, Integer *irow, Integer *icol,
                    Integer *nb, Integer *istb, Integer *indb, Integer *lindb,
                    Integer *nover, Integer *iwork);

int main(void)
{

  /* Scalars */
  double dtolg, rnorm, tol;
  Integer i, itn, k, la, lfillg, lindb, liwork, m, maxitn, mb, n, nb, nnz;
  Integer nnzc, nover, exit_status = 0;
  Nag_SparseNsym_Method method;
  Nag_SparseNsym_Piv pstrag;
  Nag_SparseNsym_Fact milug;

  /* Arrays */
  char nag_enum_arg[40];
  double *dtol;
  Complex *a, *b, *x;
  Integer *icol, *idiag, *indb, *ipivp, *ipivq, *irow, *istb, *istr, *iwork;
  Integer *lfill, *npivm;
  Nag_SparseNsym_Piv *pstrat;
  Nag_SparseNsym_Fact *milu;

  /* Nag Types */
  NagError fail;

  /* Print example header */
  printf("nag_sparse_nherm_precon_bdilu_solve (f11duc) Example Program ");
  printf("Results\n\n");

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

  /* Get the matrix order and number of nonzero entries. */
  scanf("%" NAG_IFMT " %*[^\n]", &n);
  scanf("%" NAG_IFMT " %*[^\n]", &nnz);

  la = 20 * nnz;
  lindb = 3 * n;
  liwork = 9 * n + 3;

  /* Allocate arrays */
  b = NAG_ALLOC(n, Complex);
  x = NAG_ALLOC(n, Complex);

  a = NAG_ALLOC(la, Complex);
  irow = NAG_ALLOC(la, Integer);
  icol = NAG_ALLOC(la, Integer);

  idiag = NAG_ALLOC(lindb, Integer);
  indb = NAG_ALLOC(lindb, Integer);
  ipivp = NAG_ALLOC(lindb, Integer);
  ipivq = NAG_ALLOC(lindb, Integer);
  istr = NAG_ALLOC(lindb + 1, Integer);

  iwork = NAG_ALLOC(liwork, Integer);

  if ((!b) || (!x) || (!a) || (!irow) || (!icol) || (!idiag) || (!indb) ||
      (!ipivp) || (!ipivq) || (!istr) || (!iwork)) {
    printf("Allocation failure!\n");
    exit_status = -1;
    goto END;
  }

  /* Initialize arrays */
  for (i = 0; i < n; i++) {
    b[i].re = 0.0;
    b[i].im = 0.0;
    x[i].re = 0.0;
    x[i].im = 0.0;
  }

  for (i = 0; i < la; i++) {
    a[i].re = 0.0;
    a[i].im = 0.0;
    irow[i] = 0;
    icol[i] = 0;
  }

  for (i = 0; i < lindb; i++) {
    indb[i] = 0;
    ipivp[i] = 0;
    ipivq[i] = 0;
    istr[i] = 0;
    idiag[i] = 0;
  }
  istr[lindb] = 0;

  for (i = 0; i < liwork; i++) {
    iwork[i] = 0;
  }

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

  /* Read the RHS b */
  for (i = 0; i < n; i++) {
    scanf(" (%lf, %lf)", &b[i].re, &b[i].im);
  }
  scanf("%*[^\n] ");

  /* nag_enum_name_to_value (x04nac): Converts NAG enum member name to value */
  scanf("%39s %*[^\n]", nag_enum_arg);
  method = (Nag_SparseNsym_Method) nag_enum_name_to_value(nag_enum_arg);

  /* Read algorithmic parameters */
  scanf("%" NAG_IFMT " %lf %*[^\n]", &lfillg, &dtolg);

  /* nag_enum_name_to_value (x04nac): Converts NAG enum member name to value */
  scanf("%39s %*[^\n]", nag_enum_arg);
  pstrag = (Nag_SparseNsym_Piv) nag_enum_name_to_value(nag_enum_arg);

  /* nag_enum_name_to_value (x04nac): Converts NAG enum member name to value */
  scanf("%39s %*[^\n]", nag_enum_arg);
  milug = (Nag_SparseNsym_Fact) nag_enum_name_to_value(nag_enum_arg);

  /* Read algorithmic parameters */
  scanf("%" NAG_IFMT " %lf %" NAG_IFMT " %*[^\n]", &m, &tol, &maxitn);
  scanf("%" NAG_IFMT " %" NAG_IFMT " %*[^\n]", &nb, &nover);

  if (nb < 1)
    {
      printf("Value read for nb is out of range\n");
      exit_status = -4;
      goto END;
    }

  /* Allocate arrays */
  dtol = NAG_ALLOC(nb, double);
  istb = NAG_ALLOC(nb + 1, Integer);
  lfill = NAG_ALLOC(nb, Integer);
  npivm = NAG_ALLOC(nb, Integer);

  pstrat = (Nag_SparseNsym_Piv *) NAG_ALLOC(nb, Nag_SparseNsym_Piv);
  milu = (Nag_SparseNsym_Fact *) NAG_ALLOC(nb, Nag_SparseNsym_Fact);

  if ((!dtol) || (!istb) || (!lfill) || (!npivm) || (!pstrat) || (!milu)) {
    printf("Allocation failure!\n");
    exit_status = -1;
    goto END;
  }

  /* Initialize arrays */
  for (i = 0; i < nb; i++) {
    dtol[i] = 0.0;
    istb[i] = 0;
    lfill[i] = 0;
    npivm[i] = 0;
  }
  istb[nb] = 0;

  /* Define diagonal block indices.
   * In this example use blocks of MB consecutive rows and initialize
   * assuming no overlap.
   */
  mb = (n + nb - 1) / nb;
  for (k = 0; k < nb; k++) {
    istb[k] = k * mb + 1;
  }
  istb[nb] = n + 1;

  for (i = 0; i < n; i++) {
    indb[i] = i + 1;
  }

  /* Modify INDB and ISTB to account for overlap. */
  overlap(&n, &nnz, irow, icol, &nb, istb, indb, &lindb, &nover, iwork);

  /* Set algorithmic parameters for each block from global values */
  for (k = 0; k < nb; k++) {
    lfill[k] = lfillg;
    dtol[k] = dtolg;
    pstrat[k] = pstrag;
    milu[k] = milug;
  }

  /* Initialize fail */
  INIT_FAIL(fail);

  /* Calculate factorization
   * 
   * nag_sparse_nherm_precon_bdilu (f11dtc). Calculates incomplete LU
   * factorization of local or overlapping diagonal blocks, mostly used
   * as incomplete LU preconditioner for complex sparse matrix.
   */
  nag_sparse_nherm_precon_bdilu(n, nnz, a, la, irow, icol, nb, istb, indb,
                                lindb, lfill, dtol, pstrat, milu, ipivp,
                                ipivq, istr, idiag, &nnzc, npivm, &fail);

  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_nherm_precon_bdilu (f11dtc).\n%s\n\n",
           fail.message);
    exit_status = -2;
    goto END;
  }

  /* Initialize fail */
  INIT_FAIL(fail);

  /* Solve Ax = b using nag_sparse_nherm_precon_bdilu_solve (f11duc)
   *
   * nag_sparse_nherm_precon_bdilu_solve (f11duc): Solves complex sparse
   * nonsymmetric linear system, using block-jacobi preconditioner
   * generated by f11dtc.
   */
  nag_sparse_nherm_precon_bdilu_solve(method, n, nnz, a, la, irow, icol, nb,
                                      istb, indb, lindb, ipivp, ipivq, istr,
                                      idiag, b, m, tol, maxitn, x, &rnorm,
                                      &itn, &fail);

  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_nherm_precon_bdilu_solve (f11duc).\n\n%s",
           fail.message);
    exit_status = -3;
    goto END;
  }

  /* Print output */
  printf(" Converged in %9" NAG_IFMT " iterations\n", itn);
  printf(" Final residual norm = %15.6E\n", rnorm);

  /* Output x */
  printf(" Solution vector  x\n");
  printf(" -----------------------\n");
  for (i = 0; i < n; i++) {
    printf(" ( %f,  %f )\n", x[i].re, x[i].im);
  }
  printf("\n");

END:
  NAG_FREE(b);
  NAG_FREE(x);
  NAG_FREE(a);
  NAG_FREE(irow);
  NAG_FREE(icol);
  NAG_FREE(idiag);
  NAG_FREE(indb);
  NAG_FREE(ipivp);
  NAG_FREE(ipivq);
  NAG_FREE(istr);
  NAG_FREE(iwork);
  NAG_FREE(dtol);
  NAG_FREE(istb);
  NAG_FREE(lfill);
  NAG_FREE(npivm);
  NAG_FREE(pstrat);
  NAG_FREE(milu);

  return exit_status;
}

/* ************************************************************************** */

static void overlap(Integer *n, Integer *nnz, Integer *irow, Integer *icol,
                    Integer *nb, Integer *istb, Integer *indb, Integer *lindb,
                    Integer *nover, Integer *iwork)
{
  /* Purpose
   * =======
   *
   * This routine takes a set of row indices INDB defining the diagonal blocks
   * to be used in nag_sparse_nherm_precon_bdilu (f11dtc) to define a block
   * Jacobi or additive Schwarz preconditioner, and expands them to allow for
   * NOVER levels of overlap.
   *
   * The pointer array ISTB is also updated accordingly, so that the returned
   * values of ISTB and INDB can be passed on to
   * nag_sparse_nherm_precon_bdilu (f11dtc) to define overlapping diagonal
   * blocks.
   *
   * ----------------------------------------------------------------------- */

  /* Scalars */
  Integer i, ik, ind, iover, j, k, l, n21, nadd, row;

  /* Find the number of nonzero elements in each row of the matrix A, and start
   * address of each row. Store the start addresses in iwork(n,...,2*n-1).
   */

  for (i = 0; i < (*n); i++) {
    iwork[i] = 0;
  }

  for (i = 0; i < (*nnz); i++) {
    iwork[irow[i] - 1] = iwork[irow[i] - 1] + 1;
  }
  iwork[(*n)] = 1;

  for (i = 0; i < (*n); i++) {
    iwork[(*n) + i + 1] = iwork[(*n) + i] + iwork[i];
  }

  /* Loop over blocks. */
  for (k = 0; k < (*nb); k++) {

    /* Initialize marker array. */
    for (j = 0; j < (*n); j++) {
      iwork[j] = 0;
    }

    /* Mark the rows already in block K in the workspace array. */
    for (l = istb[k]; l < istb[k + 1]; l++) {
      iwork[indb[l - 1] - 1] = 1;
    }

    /* Loop over levels of overlap. */
    for (iover = 1; iover <= (*nover); iover++) {

      /* Initialize counter of new row indices to be added. */
      ind = 0;

      /* Loop over the rows currently in the diagonal block. */
      for (l = istb[k]; l < istb[k + 1]; l++) {
        row = indb[l - 1];

        /* Loop over nonzero elements in row ROW. */
        for (i = iwork[(*n) + row - 1]; i < iwork[(*n) + row]; i++) {

          /* If the column index of the nonzero element is not in the
           * existing set for this block, store it to be added later, and
           * mark it in the marker array.
           */
          if (iwork[icol[i - 1] - 1] == 0) {
            iwork[icol[i - 1] - 1] = 1;
            iwork[2 * (*n) + 1 + ind] = icol[i - 1];
            ind = ind + 1;
          }
        }
      }

      /* Shift the indices in INDB and add the new entries for block K.
       * Change ISTB accordingly.
       */
      nadd = ind;
      if (istb[(*nb)] + nadd - 1 > (*lindb)) {
        printf("**** lindb too small, lindb = %" NAG_IFMT " ****\n", *lindb);
        exit(-1);
      }

      for (i = istb[(*nb)] - 1; i >= istb[k + 1]; i--) {
        indb[i + nadd - 1] = indb[i - 1];
      }

      n21 = 2 * (*n) + 1;
      ik = istb[k + 1] - 1;

      for (j = 0; j < nadd; j++) {
        indb[ik + j] = iwork[n21 + j];
      }

      for (j = k + 1; j < (*nb) + 1; j++) {
        istb[j] = istb[j] + nadd;
      }
    }
  }
  return;
}