1 static char help[] = "Demonstrates the use of fast Richardson for SOR. And tests\n\
2 the MatSOR() routines.\n\n";
3
4 #include <petscpc.h>
5
main(int argc,char ** args)6 int main(int argc, char **args)
7 {
8 Mat mat;
9 Vec b, u;
10 PC pc;
11 PetscInt n = 5, i, col[3];
12 PetscScalar value[3];
13
14 PetscFunctionBeginUser;
15 PetscCall(PetscInitialize(&argc, &args, NULL, help));
16 /* Create vectors */
17 PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &b));
18 PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &u));
19
20 /* Create and assemble matrix */
21 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, n, n, NULL, &mat));
22 value[0] = -1.0;
23 value[1] = 2.0;
24 value[2] = -1.0;
25 for (i = 1; i < n - 1; i++) {
26 col[0] = i - 1;
27 col[1] = i;
28 col[2] = i + 1;
29 PetscCall(MatSetValues(mat, 1, &i, 3, col, value, INSERT_VALUES));
30 }
31 i = n - 1;
32 col[0] = n - 2;
33 col[1] = n - 1;
34 PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES));
35 i = 0;
36 col[0] = 0;
37 col[1] = 1;
38 value[0] = 2.0;
39 value[1] = -1.0;
40 PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES));
41 PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
42 PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
43
44 /* Create PC context and set up data structures */
45 PetscCall(PCCreate(PETSC_COMM_WORLD, &pc));
46 PetscCall(PCSetType(pc, PCSOR));
47 PetscCall(PCSetFromOptions(pc));
48 PetscCall(PCSetOperators(pc, mat, mat));
49 PetscCall(PCSetUp(pc));
50
51 value[0] = 1.0;
52 for (i = 0; i < n; i++) {
53 PetscCall(VecSet(u, 0.0));
54 PetscCall(VecSetValues(u, 1, &i, value, INSERT_VALUES));
55 PetscCall(VecAssemblyBegin(u));
56 PetscCall(VecAssemblyEnd(u));
57 PetscCall(PCApply(pc, u, b));
58 PetscCall(VecView(b, PETSC_VIEWER_STDOUT_SELF));
59 }
60
61 /* Free data structures */
62 PetscCall(MatDestroy(&mat));
63 PetscCall(PCDestroy(&pc));
64 PetscCall(VecDestroy(&u));
65 PetscCall(VecDestroy(&b));
66 PetscCall(PetscFinalize());
67 return 0;
68 }
69
70 /*TEST
71
72 test:
73
74 TEST*/
75