/* nag_sparse_nherm_jacobi (f11dxc) Example Program.
*
* Copyright 2014 Numerical Algorithms Group.
*
* Mark 23, 2011.
*/
#include <nag.h>
#include <nag_stdlib.h>
#include <naga02.h>
#include <nagf11.h>
int main(void)
{
/* Scalars */
Integer exit_status = 0;
double anorm, sigmax, stplhs, stprhs, tol;
Integer i, irevcm, iterm, itn, lwork, lwreq, m, maxitn,
monit, n, niter, nnz;
/* Arrays */
char nag_enum_arg[100];
Complex *a = 0, *b = 0, *diag = 0, *work = 0, *x = 0;
double *wgt = 0;
Integer *icol = 0, *irow = 0;
/* NAG types */
Nag_InitializeA init;
Nag_SparseNsym_Method method;
Nag_SparseNsym_PrecType precon;
Nag_NormType norm;
Nag_SparseNsym_Weight weight;
NagError fail, fail1;
INIT_FAIL(fail);
INIT_FAIL(fail1);
printf("nag_sparse_nherm_jacobi (f11dxc) Example Program Results\n");
/* Skip heading in data file*/
scanf("%*[^\n]");
scanf("%ld%*[^\n]", &n);
scanf("%ld%*[^\n]", &nnz);
lwork = 300;
if (
!(a = NAG_ALLOC(nnz, Complex)) ||
!(b = NAG_ALLOC(n, Complex)) ||
!(diag = NAG_ALLOC(n, Complex)) ||
!(work = NAG_ALLOC(lwork, Complex)) ||
!(x = NAG_ALLOC(n, Complex)) ||
!(wgt = NAG_ALLOC(n, double)) ||
!(icol = NAG_ALLOC(nnz, Integer)) ||
!(irow = NAG_ALLOC(nnz, 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("%ld%*[^\n]", &iterm);
scanf("%ld%lf%ld%*[^\n]", &m, &tol, &maxitn);
scanf("%ld%*[^\n]", &monit);
/* Read the parameters for the preconditioner*/
scanf("%ld%*[^\n]", &niter);
anorm = 0.0;
sigmax = 0.0;
/* Read the non-zero elements of the matrix A*/
for (i = 0; i < nnz; i++)
scanf(" ( %lf , %lf ) %ld%ld%*[^\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]");
/* Call to initialize the solver */
/* nag_sparse_nherm_basic_setup (f11brc).
* Complex sparse non-Hermitian linear systems, setup
*/
nag_sparse_nherm_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_nherm_basic_setup (f11brc).\n%s\n",
fail.message);
exit_status = 1;
goto END;
}
/* Call solver repeatedly to solve the equations.
* Note: the arrays b and x are overwritten; on final exit, x will
* contain the solution and b the residual vector.
*/
irevcm = 0;
init =Nag_InitializeI;
while (irevcm != 4) {
/* nag_sparse_nherm_basic_solver (f11bsc)
* Complex sparse non-Hermitian linear systems, preconditioned RGMRES, CGS,
* Bi-CGSTAB or TFQMR method
*/
nag_sparse_nherm_basic_solver(&irevcm, x, b, wgt, work, lwreq, &fail);
switch (irevcm) {
case -1:
/* nag_sparse_nherm_matvec (f11xnc)
* Complex sparse non-Hermitian matrix vector multiply
*/
nag_sparse_nherm_matvec(Nag_ConjTrans, n, nnz, a, irow, icol,
Nag_SparseNsym_NoCheck, x, b, &fail1);
break;
case 1:
nag_sparse_nherm_matvec(Nag_NoTrans, n, nnz, a, irow, icol,
Nag_SparseNsym_NoCheck, x, b, &fail1);
break;
case 2:
/* nag_sparse_nherm_jacobi (f11dxc).
* Complex sparse nonsymmetric linear systems, line Jacobi preconditioner
*/
nag_sparse_nherm_jacobi(Nag_SparseNsym_StoreCS, Nag_NoTrans, init,
niter, n, nnz, a, irow, icol,
Nag_SparseNsym_Check, x, b, diag, &fail1);
init = Nag_InputA;
break;
case 3:
/* nag_sparse_nherm_basic_diagnostic (f11btc)
* Complex sparse nonhermitian linear systems, diagnostic
*/
nag_sparse_nherm_basic_diagnostic(&itn, &stplhs, &stprhs, &anorm,
&sigmax, work, lwreq,&fail1);
if (fail1.code == NE_NOERROR)
printf("Monitoring at iteration no.%4ld residual %14.4e\n",
itn, stplhs);
}
if (fail1.code != NE_NOERROR) irevcm = 6;
}
if (fail.code != NE_NOERROR)
{
printf("Error from nag_sparse_nherm_basic_solver (f11bsc)\n%s\n",
fail.message);
exit_status = 2;
goto END;
}
/* Obtain information about the computation using
* nag_sparse_nherm_basic_diagnostic (f11btc).
* Complex sparse Hermitian linear systems, diagnostic.
*/
nag_sparse_nherm_basic_diagnostic(&itn, &stplhs, &stprhs, &anorm, &sigmax,
work, lwreq, &fail);
if (fail.code != NE_NOERROR)
{
printf("Error from nag_sparse_nherm_basic_diagnostic (f11btc)\n%s\n",
fail.message);
exit_status = 3;
goto END;
}
/* Print the output data*/
printf("\nFinal Results\n");
printf("Number of iterations for convergence: %4ld \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%20s%29s\n","Solution vector","Residual vector");
for (i = 0; i < n; i++)
printf("(%13.4e, %13.4e) (%13.4e, %13.4e)\n", x[i].re, x[i].im, b[i].re,
b[i].im);
END:
NAG_FREE(a);
NAG_FREE(b);
NAG_FREE(diag);
NAG_FREE(work);
NAG_FREE(x);
NAG_FREE(wgt);
NAG_FREE(icol);
NAG_FREE(irow);
return exit_status;
}