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