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

NAG CL Interface Introduction
Example description
/* nag_sparse_complex_herm_solve_ilu (f11jqc) Example Program.
 *
 * Copyright 2023 Numerical Algorithms Group.
 *
 * Mark 29.1, 2023.
 */

#include <nag.h>

int main(void) {
  /* Scalars */
  Integer exit_status = 0;
  double dscale, dtol, rnorm, tol;
  Integer i, itn, la, lfill, maxitn, n, nnz, nnzc, npivm;
  /* Arrays */
  char nag_enum_arg[40];
  Complex *a = 0, *b = 0, *x = 0;
  Integer *icol = 0, *ipiv = 0, *irow = 0, *istr = 0;
  /* NAG types */
  Nag_SparseSym_Method method;
  Nag_SparseSym_Piv pstrat;
  Nag_SparseSym_Fact mic;
  NagError fail;

  INIT_FAIL(fail);

  printf(
      "nag_sparse_complex_herm_solve_ilu (f11jqc) Example Program Results\n\n");
  /* Skip heading in data file */
  scanf("%*[^\n]");
  /* Read algorithmic parameters */
  scanf("%" NAG_IFMT "%*[^\n]", &n);
  scanf("%" NAG_IFMT "%*[^\n]", &nnz);

  /* Allocate memory */
  la = 3 * nnz;
  if (!(a = NAG_ALLOC(la, Complex)) || !(b = NAG_ALLOC(n, Complex)) ||
      !(x = NAG_ALLOC(n, Complex)) || !(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;
  }
  scanf("%39s%*[^\n]", nag_enum_arg);
  /* nag_enum_name_to_value (x04nac).
   * Converts NAG enum member name to value
   */
  method = (Nag_SparseSym_Method)nag_enum_name_to_value(nag_enum_arg);
  scanf("%" NAG_IFMT "%lf%*[^\n]", &lfill, &dtol);
  scanf("%39s%*[^\n]", nag_enum_arg);
  mic = (Nag_SparseSym_Fact)nag_enum_name_to_value(nag_enum_arg);
  scanf("%lf%*[^\n]", &dscale);
  scanf("%39s%*[^\n]", nag_enum_arg);
  pstrat = (Nag_SparseSym_Piv)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 Cholesky factorization of complex sparse Hermitian
   * matrix using nag_sparse_complex_herm_precon_ichol (f11jnc).
   */
  nag_sparse_complex_herm_precon_ichol(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_complex_herm_precon_ichol (f11jnc).\n%s\n",
           fail.message);
    exit_status = 1;
    goto END;
  }
  /* Solve Linear System. */
  /* nag_sparse_complex_herm_solve_ilu (f11jqc).
   * Solution of complex sparse Hermitian linear system, conjugate
   * gradient/Lanczos method, preconditioner computed by f11jnc
   */
  nag_sparse_complex_herm_solve_ilu(method, n, nnz, a, la, irow, icol, ipiv,
                                    istr, b, tol, maxitn, x, &rnorm, &itn,
                                    &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_complex_herm_solve_ilu.(f11jqc)\n%s\n",
           fail.message);
    exit_status = 2;
    goto END;
  }
  printf("Converged in  %10" NAG_IFMT " iterations \n", itn);
  printf("Final residual norm = %10.3e\n\n", rnorm);
  printf("     Converged Solution\n");
  /* Output x */
  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(ipiv);
  NAG_FREE(irow);
  NAG_FREE(istr);
  return exit_status;
}