NAG Library Manual, Mark 28.3
Interfaces:  FL   CL   CPP   AD 

NAG CL Interface Introduction
Example description
/* nag_sparse_complex_gen_solve_ilu (f11dqc) Example Program.
 *
 * Copyright 2022 Numerical Algorithms Group.
 *
 * Mark 28.3, 2022.
 */
#include <nag.h>
int main(void) {
  /* Scalars */
  Integer exit_status = 0;
  double dtol, rnorm, tol;
  Integer i, itn, la, lfill, m, maxitn, n, nnz, nnzc, npivm;
  /* Arrays */
  Complex *a = 0, *b = 0, *x = 0;
  Integer *icol = 0, *idiag = 0, *ipivp = 0, *ipivq = 0, *irow = 0, *istr = 0;
  char nag_enum_arg[40];
  /* NAG types */
  Nag_SparseNsym_Method method;
  Nag_SparseNsym_Piv pstrat;
  Nag_SparseNsym_Fact milu;
  NagError fail;

  INIT_FAIL(fail);

  printf(
      "nag_sparse_complex_gen_solve_ilu (f11dqc) Example Program Results\n\n");

  /* Skip heading in data file */
  scanf("%*[^\n]");
  /* Read algorithmic parameters */
  scanf("%" NAG_IFMT "%" NAG_IFMT "%*[^\n]", &n, &m);
  scanf("%" NAG_IFMT "%*[^\n]", &nnz);
  la = 2 * nnz;
  if (!(a = NAG_ALLOC((la), Complex)) || !(b = NAG_ALLOC((n), Complex)) ||
      !(x = NAG_ALLOC((n), Complex)) || !(icol = NAG_ALLOC((la), Integer)) ||
      !(idiag = NAG_ALLOC((n), Integer)) ||
      !(ipivp = NAG_ALLOC((n), Integer)) ||
      !(ipivq = NAG_ALLOC((n), Integer)) ||
      !(irow = NAG_ALLOC((la), Integer)) ||
      !(istr = NAG_ALLOC((n + 1), Integer))) {
    printf("Allocation failure\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
   */
  method = (Nag_SparseNsym_Method)nag_enum_name_to_value(nag_enum_arg);
  scanf("%" NAG_IFMT "%lf%*[^\n]", &lfill, &dtol);
  scanf("%39s%*[^\n]", nag_enum_arg);
  pstrat = (Nag_SparseNsym_Piv)nag_enum_name_to_value(nag_enum_arg);
  scanf("%39s%*[^\n]", nag_enum_arg);
  milu = (Nag_SparseNsym_Fact)nag_enum_name_to_value(nag_enum_arg);
  scanf("%lf%" NAG_IFMT "%*[^\n]", &tol, &maxitn);
  /* Read 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 rhs 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 LU factorization */
  /* nag_sparse_complex_gen_precon_ilu (f11dnc)
   * Complex sparse non-Hermitian linear systems, incomplete LU factorization
   */
  nag_sparse_complex_gen_precon_ilu(n, nnz, a, la, irow, icol, lfill, dtol,
                                    pstrat, milu, ipivp, ipivq, istr, idiag,
                                    &nnzc, &npivm, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_complex_gen_precon_ilu (f11dnc).\n%s\n",
           fail.message);
    exit_status = 1;
    goto END;
  }
  /* solve ax = b */
  /* nag_sparse_complex_gen_solve_ilu (f11dqc).
   * Solution of complex sparse non-Hermitian linear system, RGMRES, CGS,
   * Bi-CGSTAB or TFQMR method, preconditioner computed by
   * nag_sparse_complex_gen_precon_ilu (f11dnc) (Black Box).
   */
  nag_sparse_complex_gen_solve_ilu(method, n, nnz, a, la, irow, icol, ipivp,
                                   ipivq, istr, idiag, b, m, tol, maxitn, x,
                                   &rnorm, &itn, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_complex_gen_solve_ilu (f11dqc).\n%s\n",
           fail.message);
    exit_status = 2;
    goto END;
  }
  printf("Converged in%12" NAG_IFMT " iterations\n", itn);
  printf("Final residual norm =%11.3e\n\n", rnorm);
  /* Output x */
  printf("%16s\n", "Solution");
  for (i = 0; i < n; i++)
    printf(" (%13.4e, %13.4e) \n", x[i].re, x[i].im);

END:
  NAG_FREE(a);
  NAG_FREE(b);
  NAG_FREE(x);
  NAG_FREE(icol);
  NAG_FREE(idiag);
  NAG_FREE(ipivp);
  NAG_FREE(ipivq);
  NAG_FREE(irow);
  NAG_FREE(istr);
  return exit_status;
}