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

NAG CL Interface Introduction
Example description
/* nag_sparse_real_gen_solve_jacssor (f11dec) Example Program.
 *
 * Copyright 2024 Numerical Algorithms Group.
 *
 * Mark 30.1, 2024.
 */

#include <nag.h>
#include <stdio.h>

int main(void) {
  double *a = 0, *b = 0, *x = 0;
  double omega;
  double rnorm;
  double tol;
  Integer exit_status = 0;
  Integer *icol = 0, *irow = 0;
  Integer i, m, n;
  Integer maxitn, itn;
  Integer nnz;
  char nag_enum_arg[40];
  Nag_SparseNsym_Method method;
  Nag_SparseNsym_PrecType precon;
  Nag_Sparse_Comm comm;
  NagError fail;

  INIT_FAIL(fail);

  printf(
      "nag_sparse_real_gen_solve_jacssor (f11dec) Example Program Results\n");
  /* Skip heading in data file */
  scanf("%*[^\n]");
  scanf("%" NAG_IFMT "%*[^\n]", &n);
  scanf("%" NAG_IFMT "%*[^\n]", &nnz);
  scanf("%39s", 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("%39s%*[^\n]", nag_enum_arg);
  precon = (Nag_SparseNsym_PrecType)nag_enum_name_to_value(nag_enum_arg);
  scanf("%lf%*[^\n]", &omega);
  scanf("%" NAG_IFMT "%lf%" NAG_IFMT "%*[^\n]", &m, &tol, &maxitn);

  x = NAG_ALLOC(n, double);
  b = NAG_ALLOC(n, double);
  a = NAG_ALLOC(nnz, double);
  irow = NAG_ALLOC(nnz, Integer);
  icol = NAG_ALLOC(nnz, Integer);
  if (!irow || !icol || !a || !x || !b) {
    printf("Allocation failure\n");
    exit_status = 1;
    goto END;
  }

  /* Read the matrix a */

  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]");

  /* Solve Ax = b using nag_sparse_real_gen_solve_jacssor (f11dec) */

  /* nag_sparse_real_gen_solve_jacssor (f11dec).
   * Solver with no Jacobi/SSOR preconditioning (nonsymmetric)
   */
  nag_sparse_real_gen_solve_jacssor(method, precon, n, nnz, a, irow, icol,
                                    omega, b, m, tol, maxitn, x, &rnorm, &itn,
                                    &comm, &fail);

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

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

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

  return exit_status;
}