xref: /petsc/src/ksp/pc/tests/ex4.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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   PetscFunctionBeginUser;
16   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
17   /* Create vectors */
18   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&b));
19   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&u));
20 
21   /* Create and assemble matrix */
22   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,n,n,NULL,&mat));
23   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
24   for (i=1; i<n-1; i++) {
25     col[0] = i-1; col[1] = i; col[2] = i+1;
26     PetscCall(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES));
27   }
28   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
29   PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
30   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
31   PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
32   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
33   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
34 
35   /* Create PC context and set up data structures */
36   PetscCall(PCCreate(PETSC_COMM_WORLD,&pc));
37   PetscCall(PCSetType(pc,PCSOR));
38   PetscCall(PCSetFromOptions(pc));
39   PetscCall(PCSetOperators(pc,mat,mat));
40   PetscCall(PCSetUp(pc));
41 
42   value[0] = 1.0;
43   for (i=0; i<n; i++) {
44     PetscCall(VecSet(u,0.0));
45     PetscCall(VecSetValues(u,1,&i,value,INSERT_VALUES));
46     PetscCall(VecAssemblyBegin(u));
47     PetscCall(VecAssemblyEnd(u));
48     PetscCall(PCApply(pc,u,b));
49     PetscCall(VecView(b,PETSC_VIEWER_STDOUT_SELF));
50   }
51 
52   /* Free data structures */
53   PetscCall(MatDestroy(&mat));
54   PetscCall(PCDestroy(&pc));
55   PetscCall(VecDestroy(&u));
56   PetscCall(VecDestroy(&b));
57   PetscCall(PetscFinalize());
58   return 0;
59 }
60 
61 /*TEST
62 
63    test:
64 
65 TEST*/
66