Example description
/* nag_sparse_nsym_jacobi (f11dkc) Example Program.
 *
 * Copyright 2017 Numerical Algorithms Group.
 *
 * Mark 26.2, 2017.
 */
#include <nag.h>
#include <nag_stdlib.h>
#include <nagf11.h>
int main(void)
{
  /* Scalars */
  Integer exit_status = 0;
  double anorm, sigmax, stplhs, stprhs, tol;
  Integer i, irevcm, iterm, itn, lwork, lwreq, m, maxitn, monit,
         n, niter, nnz;
  /* Arrays */
  char nag_enum_arg[100];
  double *a = 0, *b = 0, *diag = 0, *wgt = 0, *work = 0, *x = 0;
  Integer *icol = 0, *irow = 0;
  /* NAG types */
  Nag_InitializeA init;
  Nag_SparseNsym_Method method;
  Nag_SparseNsym_PrecType precon;
  Nag_NormType norm;
  Nag_SparseNsym_Weight weight;
  NagError fail, fail1;

  INIT_FAIL(fail);
  INIT_FAIL(fail1);

  printf("nag_sparse_nsym_jacobi (f11dkc) Example Program Results \n");
  /* Skip heading in data file */
  scanf("%*[^\n]");
  scanf("%" NAG_IFMT "%*[^\n]", &n);
  scanf("%" NAG_IFMT "%*[^\n]", &nnz);
  lwork = 200;
  if (!(a = NAG_ALLOC((nnz), double)) ||
      !(b = NAG_ALLOC((n), double)) ||
      !(diag = NAG_ALLOC((n), double)) ||
      !(wgt = NAG_ALLOC((n), double)) ||
      !(work = NAG_ALLOC((lwork), double)) ||
      !(x = NAG_ALLOC((n), double)) ||
      !(icol = NAG_ALLOC((nnz), Integer)) ||
      !(irow = NAG_ALLOC((nnz), Integer))
         )
  {
    printf("Allocation failure\n");
    exit_status = -1;
    goto END;
  }
  /* Read or initialize the parameters for the iterative solver */
  scanf("%99s%*[^\n]", nag_enum_arg);
  /* nag_enum_name_to_value (x04nac).
   * Converts NAG enum member name to value
   */
  method = (Nag_SparseNsym_Method) nag_enum_name_to_value(nag_enum_arg);
  scanf("%99s%*[^\n]", nag_enum_arg);
  precon = (Nag_SparseNsym_PrecType) 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_SparseNsym_Weight) nag_enum_name_to_value(nag_enum_arg);
  scanf("%" NAG_IFMT "%*[^\n]", &iterm);
  scanf("%" NAG_IFMT "%lf%" NAG_IFMT "%*[^\n]", &m, &tol, &maxitn);
  scanf("%" NAG_IFMT "%*[^\n]", &monit);

  /* Read the parameters for the preconditioner */
  scanf("%" NAG_IFMT "%*[^\n]", &niter);
  anorm = 0.0;
  sigmax = 0.0;

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

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

  /* nag_sparse_nsym_basic_setup (f11bdc)
   * Real sparse nonsymmetric linear systems, setup routine
   */
  nag_sparse_nsym_basic_setup(method, precon, norm, weight, iterm, n, m, tol,
                              maxitn, anorm, sigmax, monit, &lwreq, work,
                              lwork, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_nsym_basic_setup (f11bdc).\n%s\n",
           fail.message);
    exit_status = 1;
    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;
  init = Nag_InitializeI;
  while (irevcm != 4) {
    /* nag_sparse_nsym_basic_solver (f11bec)
     * Real sparse nonsymmetric linear systems, solver routine
     * preconditioned RGMRES, CGS, Bi-CGSTAB or TFQMR method
     */
    nag_sparse_nsym_basic_solver(&irevcm, x, b, wgt, work, lwreq, &fail);
    switch (irevcm) {
    case -1:
      /* nag_sparse_nsym_matvec (f11xac)
       * Real sparse nonsymmetric matrix vector multiply
       */
      nag_sparse_nsym_matvec(Nag_Trans, n, nnz, a, irow, icol,
                             Nag_SparseNsym_NoCheck, x, b, &fail1);
      break;
    case 1:
      nag_sparse_nsym_matvec(Nag_NoTrans, n, nnz, a, irow, icol,
                             Nag_SparseNsym_NoCheck, x, b, &fail1);
      break;
    case 2:
      /* nag_sparse_nsym_jacobi (f11dkc)
       * Real sparse nonsymmetric linear systems, line Jacobi preconditioner
       */
      nag_sparse_nsym_jacobi(Nag_SparseNsym_StoreCS, Nag_NoTrans, init,
                             niter, n, nnz, a, irow, icol,
                             Nag_SparseNsym_Check, x, b, diag, &fail1);
      init = Nag_InputA;
      break;
    case 3:
      /* nag_sparse_nsym_basic_diagnostic (f11bfc)
       * Real sparse nonsymmetric linear systems, diagnostic for f11bec
       */
      nag_sparse_nsym_basic_diagnostic(&itn, &stplhs, &stprhs, &anorm,
                                       &sigmax, work, lwreq, &fail1);
      printf("%" NAG_IFMT " %f \n", itn, stplhs);
    }
    if (fail1.code != NE_NOERROR)
      irevcm = 6;
  }
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_nsym_basic_solver (f11bec)\n%s\n",
           fail.message);
    exit_status = 2;
    goto END;
  }

  /* Obtain information about the computation */
  nag_sparse_nsym_basic_diagnostic(&itn, &stplhs, &stprhs, &anorm, &sigmax,
                                   work, lwreq, &fail);
  /* Print the output data */
  printf("\nFinal 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\n", anorm);

  /* Output x */
  printf("%16s%16s\n", "Solution", "Residuals");
  for (i = 0; i < n; i++)
    printf("%16.4f%16.4e\n", x[i], b[i]);

END:
  NAG_FREE(a);
  NAG_FREE(b);
  NAG_FREE(diag);
  NAG_FREE(wgt);
  NAG_FREE(work);
  NAG_FREE(x);
  NAG_FREE(icol);
  NAG_FREE(irow);
  return exit_status;
}