xref: /petsc/src/ksp/pc/tests/ex3.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 
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";
5 
6 #include <petscksp.h>
7 #include <petscpc.h>
8 
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;
19 
20   PetscFunctionBeginUser;
21   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
22   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
23 
24   /* Create and initialize vectors */
25   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&b));
26   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&ustar));
27   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&u));
28   PetscCall(VecSet(ustar,1.0));
29   PetscCall(VecSet(u,0.0));
30 
31   /* Create and assemble matrix */
32   PetscCall(MatCreate(PETSC_COMM_SELF,&mat));
33   PetscCall(MatSetType(mat,MATSEQAIJ));
34   PetscCall(MatSetSizes(mat,n,n,n,n));
35   PetscCall(MatSetFromOptions(mat));
36   PetscCall(MatSeqAIJSetPreallocation(mat,3,NULL));
37   PetscCall(MatSeqBAIJSetPreallocation(mat,1,3,NULL));
38   PetscCall(MatSeqSBAIJSetPreallocation(mat,1,3,NULL));
39   PetscCall(MatSeqSELLSetPreallocation(mat,3,NULL));
40   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
41   for (i=1; i<n-1; i++) {
42     col[0] = i-1; col[1] = i; col[2] = i+1;
43     PetscCall(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES));
44   }
45   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
46   PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
47   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
48   PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
49   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
50   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
51 
52   /* Compute right-hand-side vector */
53   PetscCall(MatMult(mat,ustar,b));
54 
55   /* Create PC context and set up data structures */
56   PetscCall(PCCreate(PETSC_COMM_WORLD,&pc));
57   PetscCall(PCSetType(pc,PCNONE));
58   PetscCall(PCSetFromOptions(pc));
59   PetscCall(PCSetOperators(pc,mat,mat));
60   PetscCall(PCSetUp(pc));
61 
62   /* Create KSP context and set up data structures */
63   PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp));
64   PetscCall(KSPSetType(ksp,KSPRICHARDSON));
65   PetscCall(KSPSetFromOptions(ksp));
66   PetscCall(PCSetOperators(pc,mat,mat));
67   PetscCall(KSPSetPC(ksp,pc));
68   PetscCall(KSPSetUp(ksp));
69 
70   /* Solve the problem */
71   PetscCall(KSPGetType(ksp,&kspname));
72   PetscCall(PCGetType(pc,&pcname));
73   PetscCall(PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname));
74   PetscCall(KSPSolve(ksp,b,u));
75   PetscCall(KSPGetIterationNumber(ksp,&its));
76   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of iterations %" PetscInt_FMT "\n",its));
77 
78   /* Free data structures */
79   PetscCall(KSPDestroy(&ksp));
80   PetscCall(VecDestroy(&u));
81   PetscCall(VecDestroy(&ustar));
82   PetscCall(VecDestroy(&b));
83   PetscCall(MatDestroy(&mat));
84   PetscCall(PCDestroy(&pc));
85   PetscCall(PetscFinalize());
86   return 0;
87 }
88 
89 /*TEST
90 
91    testset:
92      args: -ksp_type gmres -ksp_monitor_short -pc_type sor -pc_sor_symmetric
93      output_file: output/ex3_1.out
94      test:
95        suffix: sor_aij
96      test:
97        suffix: sor_seqbaij
98        args: -mat_type seqbaij
99      test:
100        suffix: sor_seqsbaij
101        args: -mat_type seqbaij
102      test:
103        suffix: sor_seqsell
104        args: -mat_type seqsell
105 
106 TEST*/
107