/* nag_sparse_complex_gen_basic_setup (f11brc) Example Program.
*
* Copyright 2024 Numerical Algorithms Group.
*
* Mark 30.0, 2024.
*/
#include <nag.h>
int main(void) {
/* Scalars */
Integer exit_status = 0;
double anorm, dtol, sigmax, stplhs, stprhs, tol;
Integer i, irevcm, iterm, itn, la, lfill, lwork, lwreq, m, maxitn, monit, n,
nnz, nnzc, npivm;
/* Arrays */
char nag_enum_arg[100];
Complex *a = 0, *b = 0, *work = 0, *x = 0;
double *wgt = 0;
Integer *icol = 0, *idiag = 0, *ipivp = 0, *ipivq = 0, *irow = 0, *istr = 0;
/* NAG types */
Nag_SparseNsym_Piv pstrat;
Nag_SparseNsym_Fact milu;
Nag_SparseNsym_Method method;
Nag_SparseNsym_PrecType precon;
Nag_NormType norm;
Nag_SparseNsym_Weight weight;
Nag_TransType trans;
Nag_SparseNsym_CheckData check = Nag_SparseNsym_NoCheck;
NagError fail, fail1;
INIT_FAIL(fail);
INIT_FAIL(fail1);
printf(
"nag_sparse_complex_gen_basic_setup (f11brc) Example Program Results\n");
/* Skip heading in data file */
scanf("%*[^\n]");
scanf("%" NAG_IFMT "%*[^\n]", &n);
scanf("%" NAG_IFMT "%*[^\n]", &nnz);
la = 2 * nnz;
lwork = 200;
if (!(a = NAG_ALLOC((la), Complex)) || !(b = NAG_ALLOC((n), Complex)) ||
!(work = NAG_ALLOC((lwork), Complex)) || !(x = NAG_ALLOC((n), Complex)) ||
!(wgt = NAG_ALLOC((n), double)) || !(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;
}
/* 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 "%lf%*[^\n]", &lfill, &dtol);
scanf("%99s%*[^\n]", nag_enum_arg);
milu = (Nag_SparseNsym_Fact)nag_enum_name_to_value(nag_enum_arg);
scanf("%99s%*[^\n]", nag_enum_arg);
pstrat = (Nag_SparseNsym_Piv)nag_enum_name_to_value(nag_enum_arg);
/* Read the nonzero elements of 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 right-hand side vector B and initial approximate solution */
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);
scanf("%*[^\n]");
/* nag_sparse_complex_gen_precon_ilu (f11dnc)
* Incomplete LU factorization (non-hermitian)
*/
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;
}
/* Initialize the solver using nag_sparse_complex_gen_basic_setup (f11brc)
* Complex sparse non-Hermitian linear systems, setup routine
*/
anorm = 0.0;
sigmax = 0.0;
nag_sparse_complex_gen_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_complex_gen_basic_setup (f11brc)\n%s\n",
fail.message);
exit_status = 2;
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;
/* First call to nag_sparse_complex_gen_basic_solver (f11bsc)
* Complex sparse non-Hermitian linear systems, solver routine
* preconditioned RGMRES, CGS, Bi-CGSTAB or TFQMR method
*/
nag_sparse_complex_gen_basic_solver(&irevcm, x, b, wgt, work, lwreq, &fail);
while (irevcm != 4) {
switch (irevcm) {
case -1:
/* nag_sparse_complex_gen_matvec (f11xnc)
* Complex sparse non-Hermitian matrix vector multiply
*/
trans = Nag_ConjTrans;
nag_sparse_complex_gen_matvec(trans, n, nnz, a, irow, icol, check, x, b,
&fail1);
break;
case 1:
trans = Nag_NoTrans;
nag_sparse_complex_gen_matvec(trans, n, nnz, a, irow, icol, check, x, b,
&fail1);
break;
case 2:
/* nag_sparse_complex_gen_precon_ilu_solve (f11dpc)
* Solution of complex linear system involving incomplete LU
* preconditioning matrix
*/
trans = Nag_NoTrans;
nag_sparse_complex_gen_precon_ilu_solve(trans, n, a, la, irow, icol,
ipivp, ipivq, istr, idiag, check,
x, b, &fail1);
break;
case 3:
/* nag_sparse_complex_gen_basic_diag (f11btc)
* Complex sparse non-Hermitian linear systems, diagnostic routine
*/
nag_sparse_complex_gen_basic_diag(&itn, &stplhs, &stprhs, &anorm, &sigmax,
work, lwreq, &fail1);
printf("\nMonitoring at iteration no.%4" NAG_IFMT "\n", itn);
printf("residual norm%14.4e\n\n", stplhs);
printf(" Current Solution vector\n");
for (i = 0; i < n; i++)
printf(" (%13.4e, %13.4e)\n", x[i].re, x[i].im);
printf("\n Current Residual vector\n");
for (i = 0; i < n; i++)
printf(" (%13.4e, %13.4e)\n", b[i].re, b[i].im);
printf("\n");
}
if (fail1.code != NE_NOERROR)
irevcm = 6;
/* Next call to nag_sparse_complex_gen_basic_solver (f11bsc) */
nag_sparse_complex_gen_basic_solver(&irevcm, x, b, wgt, work, lwreq, &fail);
}
if (fail.code != NE_NOERROR) {
printf("Error from nag_sparse_complex_gen_basic_solver (f11bsc)\n%s\n",
fail.message);
exit_status = 3;
goto END;
}
/* Obtain information about the computation using
* nag_sparse_complex_gen_basic_diag (f11btc).
*/
nag_sparse_complex_gen_basic_diag(&itn, &stplhs, &stprhs, &anorm, &sigmax,
work, lwreq, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_sparse_complex_gen_basic_diag (f11btc) \n%s\n",
fail.message);
exit_status = 4;
goto END;
}
/* Print the output data */
printf("Final 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", anorm);
/* Output x */
printf("\n Solution vector\n");
for (i = 0; i < n; i++)
printf(" (%13.4e, %13.4e)\n", x[i].re, x[i].im);
printf("\n Residual vector\n");
for (i = 0; i < n; i++)
printf(" (%13.4e, %13.4e)\n", b[i].re, b[i].im);
printf("\n");
END:
NAG_FREE(a);
NAG_FREE(b);
NAG_FREE(work);
NAG_FREE(x);
NAG_FREE(wgt);
NAG_FREE(icol);
NAG_FREE(idiag);
NAG_FREE(ipivp);
NAG_FREE(ipivq);
NAG_FREE(irow);
NAG_FREE(istr);
return exit_status;
}