Actual source code: ex3.c
2: static char help[] = "Demonstrates the use of fast Richardson for SOR. And\n\
3: also tests the MatSOR() routines. Input parameters are:\n\
4: -n <n> : problem dimension\n\n";
6: #include <petscksp.h>
7: #include <petscpc.h>
9: int main(int argc, char **args)
10: {
11: Mat mat; /* matrix */
12: Vec b, ustar, u; /* vectors (RHS, exact solution, approx solution) */
13: PC pc; /* PC context */
14: KSP ksp; /* KSP context */
15: PetscInt n = 10, i, its, col[3];
16: PetscScalar value[3];
17: KSPType kspname;
18: PCType pcname;
21: PetscInitialize(&argc, &args, (char *)0, help);
22: PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL);
24: /* Create and initialize vectors */
25: VecCreateSeq(PETSC_COMM_SELF, n, &b);
26: VecCreateSeq(PETSC_COMM_SELF, n, &ustar);
27: VecCreateSeq(PETSC_COMM_SELF, n, &u);
28: VecSet(ustar, 1.0);
29: VecSet(u, 0.0);
31: /* Create and assemble matrix */
32: MatCreate(PETSC_COMM_SELF, &mat);
33: MatSetType(mat, MATSEQAIJ);
34: MatSetSizes(mat, n, n, n, n);
35: MatSetFromOptions(mat);
36: MatSeqAIJSetPreallocation(mat, 3, NULL);
37: MatSeqBAIJSetPreallocation(mat, 1, 3, NULL);
38: MatSeqSBAIJSetPreallocation(mat, 1, 3, NULL);
39: MatSeqSELLSetPreallocation(mat, 3, NULL);
40: value[0] = -1.0;
41: value[1] = 2.0;
42: value[2] = -1.0;
43: for (i = 1; i < n - 1; i++) {
44: col[0] = i - 1;
45: col[1] = i;
46: col[2] = i + 1;
47: MatSetValues(mat, 1, &i, 3, col, value, INSERT_VALUES);
48: }
49: i = n - 1;
50: col[0] = n - 2;
51: col[1] = n - 1;
52: MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES);
53: i = 0;
54: col[0] = 0;
55: col[1] = 1;
56: value[0] = 2.0;
57: value[1] = -1.0;
58: MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES);
59: MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);
60: MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);
62: /* Compute right-hand-side vector */
63: MatMult(mat, ustar, b);
65: /* Create PC context and set up data structures */
66: PCCreate(PETSC_COMM_WORLD, &pc);
67: PCSetType(pc, PCNONE);
68: PCSetFromOptions(pc);
69: PCSetOperators(pc, mat, mat);
70: PCSetUp(pc);
72: /* Create KSP context and set up data structures */
73: KSPCreate(PETSC_COMM_WORLD, &ksp);
74: KSPSetType(ksp, KSPRICHARDSON);
75: KSPSetFromOptions(ksp);
76: PCSetOperators(pc, mat, mat);
77: KSPSetPC(ksp, pc);
78: KSPSetUp(ksp);
80: /* Solve the problem */
81: KSPGetType(ksp, &kspname);
82: PCGetType(pc, &pcname);
83: PetscPrintf(PETSC_COMM_SELF, "Running %s with %s preconditioning\n", kspname, pcname);
84: KSPSolve(ksp, b, u);
85: KSPGetIterationNumber(ksp, &its);
86: PetscPrintf(PETSC_COMM_WORLD, "Number of iterations %" PetscInt_FMT "\n", its);
88: /* Free data structures */
89: KSPDestroy(&ksp);
90: VecDestroy(&u);
91: VecDestroy(&ustar);
92: VecDestroy(&b);
93: MatDestroy(&mat);
94: PCDestroy(&pc);
95: PetscFinalize();
96: return 0;
97: }
99: /*TEST
101: testset:
102: args: -ksp_type gmres -ksp_monitor_short -pc_type sor -pc_sor_symmetric
103: output_file: output/ex3_1.out
104: test:
105: suffix: sor_aij
106: test:
107: suffix: sor_seqbaij
108: args: -mat_type seqbaij
109: test:
110: suffix: sor_seqsbaij
111: args: -mat_type seqbaij
112: test:
113: suffix: sor_seqsell
114: args: -mat_type seqsell
116: TEST*/