xref: /petsc/src/ksp/pc/tests/ex2.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests PC and KSP on a tridiagonal matrix.  Note that most\n\
3c4762a1bSJed Brown users should employ the KSP interface instead of using PC directly.\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscksp.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            mat;          /* matrix */
10c4762a1bSJed Brown   Vec            b,ustar,u;  /* vectors (RHS, exact solution, approx solution) */
11c4762a1bSJed Brown   PC             pc;           /* PC context */
12c4762a1bSJed Brown   KSP            ksp;          /* KSP context */
13c4762a1bSJed Brown   PetscInt       n = 10,i,its,col[3];
14c4762a1bSJed Brown   PetscScalar    value[3];
15c4762a1bSJed Brown   PCType         pcname;
16c4762a1bSJed Brown   KSPType        kspname;
17c4762a1bSJed Brown   PetscReal      norm,tol=1000.*PETSC_MACHINE_EPSILON;
18c4762a1bSJed Brown 
19*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
20c4762a1bSJed Brown   /* Create and initialize vectors */
215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&b));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&ustar));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&u));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(ustar,1.0));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u,0.0));
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /* Create and assemble matrix */
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&mat));
29c4762a1bSJed Brown   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
30c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
31c4762a1bSJed Brown     col[0] = i-1; col[1] = i; col[2] = i+1;
325f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES));
33c4762a1bSJed Brown   }
34c4762a1bSJed Brown   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
36c4762a1bSJed Brown   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Compute right-hand-side vector */
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(mat,ustar,b));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Create PC context and set up data structures */
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PCCreate(PETSC_COMM_WORLD,&pc));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCNONE));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetFromOptions(pc));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetOperators(pc,mat,mat));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetUp(pc));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Create KSP context and set up data structures */
525f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(PETSC_COMM_WORLD,&ksp));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetType(ksp,KSPRICHARDSON));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(ksp));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetOperators(pc,mat,mat));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetPC(ksp,pc));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetUp(ksp));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Solve the problem */
605f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetType(ksp,&kspname));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PCGetType(pc,&pcname));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolve(ksp,b,u));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(u,-1.0,ustar));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(u,NORM_2,&norm));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(ksp,&its));
67c4762a1bSJed Brown   if (norm > tol) {
685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"2 norm of error %g Number of iterations %D\n",(double)norm,its));
69c4762a1bSJed Brown   }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Free data structures */
725f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&ksp));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ustar));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(PCDestroy(&pc));
78c4762a1bSJed Brown 
79*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
80*b122ec5aSJacob Faibussowitsch   return 0;
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown /*TEST
84c4762a1bSJed Brown 
85c4762a1bSJed Brown    test:
86c4762a1bSJed Brown       args: -ksp_type cg -ksp_monitor_short
87c4762a1bSJed Brown 
88c4762a1bSJed Brown TEST*/
89