Example description
/* nag_sparse_sym_chol_sol (f11jcc) Example Program.
 *
 * Copyright 2017 Numerical Algorithms Group.
 *
 * Mark 26.2, 2017.
 *
 */

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

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

  INIT_FAIL(fail);

  printf("nag_sparse_sym_chol_sol (f11jcc) Example Program Results\n");

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

  /* Read algorithmic parameters */
  scanf("%" NAG_IFMT "%*[^\n]", &n);
  scanf("%" NAG_IFMT "%*[^\n]", &nnz);
  scanf("%" NAG_IFMT "%lf%*[^\n]", &lfill, &dtol);
  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("%39s%lf%*[^\n]", nag_enum_arg, &dscale);
  mic = (Nag_SparseSym_Fact) nag_enum_name_to_value(nag_enum_arg);
  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 */

  /* Allocate memory */
  num = 2 * nnz;
  irow = NAG_ALLOC(num, Integer);
  icol = NAG_ALLOC(num, Integer);
  a = NAG_ALLOC(num, double);
  b = NAG_ALLOC(n, double);
  x = NAG_ALLOC(n, double);
  istr = NAG_ALLOC(n + 1, Integer);
  ipiv = NAG_ALLOC(num, Integer);

  if (!irow || !icol || !a || !x || !istr || !ipiv) {
    printf("Allocation failure\n");
    return EXIT_FAILURE;
  }

  for (i = 1; i <= nnz; ++i)
    scanf("%lf%" NAG_IFMT "%" NAG_IFMT "%*[^\n]", &a[i - 1], &irow[i - 1],
          &icol[i - 1]);

  /* Read right-hand side vector b and initial approximate solution x */

  for (i = 1; i <= n; ++i)
    scanf("%lf", &b[i - 1]);
  scanf(" %*[^\n]");

  for (i = 1; i <= n; ++i)
    scanf("%lf", &x[i - 1]);
  scanf("%*[^\n]");

  /* Calculate incomplete Cholesky factorization */

  /* nag_sparse_sym_chol_fac (f11jac).
   * Incomplete Cholesky factorization (symmetric)
   */
  nag_sparse_sym_chol_fac(n, nnz, &a, &num, &irow, &icol, lfill, dtol, mic,
                          dscale, pstrat, ipiv, istr, &nnzc, &npivm, &comm,
                          &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_sym_chol_fac (f11jac).\n%s\n",
           fail.message);
    exit_status = 1;
    goto END;
  }

  /* Solve Ax = b */

  /* nag_sparse_sym_chol_sol (f11jcc).
   * Solver with incomplete Cholesky preconditioning
   * (symmetric)
   */
  nag_sparse_sym_chol_sol(method, n, nnz, a, num, irow, icol, ipiv, istr, b,
                          tol, maxitn, x, &rnorm, &itn, &comm, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_sparse_sym_chol_sol (f11jcc).\n%s\n",
           fail.message);
    exit_status = 1;
    goto END;
  }

  printf(" %s%10" NAG_IFMT "%s\n", "Converged in", itn, " iterations");
  printf(" %s%16.3e\n", "Final residual norm =", rnorm);

  /* Output x */

  for (i = 1; i <= n; ++i)
    printf(" %16.4e\n", x[i - 1]);

END:
  NAG_FREE(irow);
  NAG_FREE(icol);
  NAG_FREE(a);
  NAG_FREE(b);
  NAG_FREE(x);
  NAG_FREE(ipiv);
  NAG_FREE(istr);

  return exit_status;
}