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