int SUNPCGSetPrecType(SUNLinearSolver S, int pretype) { return(SUNLinSol_PCGSetPrecType(S, pretype)); }
/* ---------------------------------------------------------------------- * SUNOCG Linear Solver Testing Routine * * We run multiple tests to exercise this solver: * 1. simple tridiagonal system (no preconditioning) * 2. simple tridiagonal system (Jacobi preconditioning) * 3. tridiagonal system w/ scale vector s (no preconditioning) * 4. tridiagonal system w/ scale vector s (Jacobi preconditioning) * * Note: We construct a tridiagonal matrix Ahat, a random solution * xhat, and a corresponding rhs vector bhat = Ahat*xhat, such * that each of these is unit-less. To test scaling, we use * the matrix * A = (S-inverse) Ahat (S-inverse), * solution vector * x = S xhat; * and construct b = A*x. Hence the linear system has both rows * and columns scaled by (S-inverse), where S is the diagonal * matrix with entries from the vector s, the 'scaling' vector * supplied to PCG having strictly positive entries. * * When this is combined with preconditioning, we construct * P \approx (A-inverse) by taking a unit-less preconditioner * Phat \approx (Ahat-inverse), and constructing the operator * P via * P = S Phat S \approx S (Ahat-inverse) S = A-inverse * We apply this via the steps: * z = Pr = S Phat S r * Since both S and Phat are diagonal matrices, this is * equivalent to * z(i) = s(i)^2 Phat(i) r(i) * --------------------------------------------------------------------*/ int main(int argc, char *argv[]) { int fails=0; /* counter for test failures */ int passfail=0; /* overall pass/fail flag */ SUNLinearSolver LS; /* linear solver object */ N_Vector xhat, x, b; /* test vectors */ UserData ProbData; /* problem data structure */ int maxl, print_timing; sunindextype i; realtype *vecdata; double tol; /* check inputs: local problem size, timing flag */ if (argc < 5) { printf("ERROR: FOUR (4) Inputs required:\n"); printf(" Problem size should be >0\n"); printf(" Maximum Krylov subspace dimension should be >0\n"); printf(" Solver tolerance should be >0\n"); printf(" timing output flag should be 0 or 1 \n"); return 1; } ProbData.N = atol(argv[1]); problem_size = ProbData.N; if (ProbData.N <= 0) { printf("ERROR: Problem size must be a positive integer\n"); return 1; } maxl = atoi(argv[2]); if (maxl <= 0) { printf("ERROR: Maximum Krylov subspace dimension must be a positive integer\n"); return 1; } tol = atof(argv[3]); if (tol <= ZERO) { printf("ERROR: Solver tolerance must be a positive real number\n"); return 1; } print_timing = atoi(argv[4]); SetTiming(print_timing); printf("\nPCG linear solver test:\n"); printf(" Problem size = %ld\n", (long int) ProbData.N); printf(" Maximum Krylov subspace dimension = %i\n", maxl); printf(" Solver Tolerance = %"GSYM"\n", tol); printf(" timing output flag = %i\n\n", print_timing); /* Create vectors */ x = N_VNew_Serial(ProbData.N); if (check_flag(x, "N_VNew_Serial", 0)) return 1; xhat = N_VNew_Serial(ProbData.N); if (check_flag(xhat, "N_VNew_Serial", 0)) return 1; b = N_VNew_Serial(ProbData.N); if (check_flag(b, "N_VNew_Serial", 0)) return 1; ProbData.d = N_VNew_Serial(ProbData.N); if (check_flag(ProbData.d, "N_VNew_Serial", 0)) return 1; ProbData.s = N_VNew_Serial(ProbData.N); if (check_flag(ProbData.s, "N_VNew_Serial", 0)) return 1; /* Fill xhat vector with uniform random data in [1,2] */ vecdata = N_VGetArrayPointer(xhat); for (i=0; i<ProbData.N; i++) vecdata[i] = ONE + urand(); /* Fill Jacobi vector with matrix diagonal */ N_VConst(FIVE, ProbData.d); /* Create PCG linear solver */ LS = SUNLinSol_PCG(x, PREC_RIGHT, maxl); fails += Test_SUNLinSolGetType(LS, SUNLINEARSOLVER_ITERATIVE, 0); fails += Test_SUNLinSolSetATimes(LS, &ProbData, ATimes, 0); fails += Test_SUNLinSolSetPreconditioner(LS, &ProbData, PSetup, PSolve, 0); fails += Test_SUNLinSolSetScalingVectors(LS, ProbData.s, NULL, 0); fails += Test_SUNLinSolInitialize(LS, 0); fails += Test_SUNLinSolSpace(LS, 0); if (fails) { printf("FAIL: SUNLinSol_PCG module failed %i initialization tests\n\n", fails); return 1; } else { printf("SUCCESS: SUNLinSol_PCG module passed all initialization tests\n\n"); } /*** Test 1: simple Poisson-like solve (no preconditioning) ***/ /* set scaling vector */ N_VConst(ONE, ProbData.s); /* Fill x vector with scaled version */ N_VProd(xhat, ProbData.s, x); /* Fill b vector with result of matrix-vector product */ fails = ATimes(&ProbData, x, b); if (check_flag(&fails, "ATimes", 1)) return 1; /* Run test with this setup */ fails += SUNLinSol_PCGSetPrecType(LS, PREC_NONE); fails += Test_SUNLinSolSetup(LS, NULL, 0); fails += Test_SUNLinSolSolve(LS, NULL, x, b, tol, 0); fails += Test_SUNLinSolLastFlag(LS, 0); fails += Test_SUNLinSolNumIters(LS, 0); fails += Test_SUNLinSolResNorm(LS, 0); fails += Test_SUNLinSolResid(LS, 0); /* Print result */ if (fails) { printf("FAIL: SUNLinSol_PCG module, problem 1, failed %i tests\n\n", fails); passfail += 1; } else { printf("SUCCESS: SUNLinSol_PCG module, problem 1, passed all tests\n\n"); } /*** Test 2: simple Poisson-like solve (Jacobi preconditioning) ***/ /* set scaling vector */ N_VConst(ONE, ProbData.s); /* Fill x vector with scaled version */ N_VProd(xhat, ProbData.s, x); /* Fill b vector with result of matrix-vector product */ fails = ATimes(&ProbData, x, b); if (check_flag(&fails, "ATimes", 1)) return 1; /* Run tests with this setup */ fails += SUNLinSol_PCGSetPrecType(LS, PREC_RIGHT); fails += Test_SUNLinSolSetup(LS, NULL, 0); fails += Test_SUNLinSolSolve(LS, NULL, x, b, tol, 0); fails += Test_SUNLinSolLastFlag(LS, 0); fails += Test_SUNLinSolNumIters(LS, 0); fails += Test_SUNLinSolResNorm(LS, 0); fails += Test_SUNLinSolResid(LS, 0); /* Print result */ if (fails) { printf("FAIL: SUNLinSol_PCG module, problem 2, failed %i tests\n\n", fails); passfail += 1; } else { printf("SUCCESS: SUNLinSol_PCG module, problem 2, passed all tests\n\n"); } /*** Test 3: Poisson-like solve w/ scaling (no preconditioning) ***/ /* set scaling vector */ vecdata = N_VGetArrayPointer(ProbData.s); for (i=0; i<ProbData.N; i++) vecdata[i] = ONE + THOUSAND*urand(); /* Fill x vector with scaled version */ N_VProd(xhat, ProbData.s, x); /* Fill b vector with result of matrix-vector product */ fails = ATimes(&ProbData, x, b); if (check_flag(&fails, "ATimes", 1)) return 1; /* Run tests with this setup */ fails += SUNLinSol_PCGSetPrecType(LS, PREC_NONE); fails += Test_SUNLinSolSetup(LS, NULL, 0); fails += Test_SUNLinSolSolve(LS, NULL, x, b, tol, 0); fails += Test_SUNLinSolLastFlag(LS, 0); fails += Test_SUNLinSolNumIters(LS, 0); fails += Test_SUNLinSolResNorm(LS, 0); fails += Test_SUNLinSolResid(LS, 0); /* Print result */ if (fails) { printf("FAIL: SUNLinSol_PCG module, problem 3, failed %i tests\n\n", fails); passfail += 1; } else { printf("SUCCESS: SUNLinSol_PCG module, problem 3, passed all tests\n\n"); } /*** Test 4: Poisson-like solve w/ scaling (Jacobi preconditioning) ***/ /* set scaling vectors */ vecdata = N_VGetArrayPointer(ProbData.s); for (i=0; i<ProbData.N; i++) vecdata[i] = ONE + THOUSAND*urand(); /* Fill x vector with scaled version */ N_VProd(xhat, ProbData.s, x); /* Fill b vector with result of matrix-vector product */ fails = ATimes(&ProbData, x, b); if (check_flag(&fails, "ATimes", 1)) return 1; /* Run tests with this setup */ fails += SUNLinSol_PCGSetPrecType(LS, PREC_RIGHT); fails += Test_SUNLinSolSetup(LS, NULL, 0); fails += Test_SUNLinSolSolve(LS, NULL, x, b, tol, 0); fails += Test_SUNLinSolLastFlag(LS, 0); fails += Test_SUNLinSolNumIters(LS, 0); fails += Test_SUNLinSolResNorm(LS, 0); fails += Test_SUNLinSolResid(LS, 0); /* Print result */ if (fails) { printf("FAIL: SUNLinSol_PCG module, problem 4, failed %i tests\n\n", fails); passfail += 1; } else { printf("SUCCESS: SUNLinSol_PCG module, problem 4, passed all tests\n\n"); } /* Free solver and vectors */ SUNLinSolFree(LS); N_VDestroy(x); N_VDestroy(xhat); N_VDestroy(b); N_VDestroy(ProbData.d); N_VDestroy(ProbData.s); return(passfail); }