xref: /petsc/src/ksp/pc/tests/ex2.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
14c4762a1bSJed Brown   PetscInt       n = 10,i,its,col[3];
15c4762a1bSJed Brown   PetscScalar    value[3];
16c4762a1bSJed Brown   PCType         pcname;
17c4762a1bSJed Brown   KSPType        kspname;
18c4762a1bSJed Brown   PetscReal      norm,tol=1000.*PETSC_MACHINE_EPSILON;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
21c4762a1bSJed Brown   /* Create and initialize vectors */
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&b));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&ustar));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&u));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(ustar,1.0));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u,0.0));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Create and assemble matrix */
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&mat));
30c4762a1bSJed Brown   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
31c4762a1bSJed Brown   for (i=1; i<n-1; i++) {
32c4762a1bSJed Brown     col[0] = i-1; col[1] = i; col[2] = i+1;
33*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES));
34c4762a1bSJed Brown   }
35c4762a1bSJed Brown   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
37c4762a1bSJed Brown   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* Compute right-hand-side vector */
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(mat,ustar,b));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Create PC context and set up data structures */
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCCreate(PETSC_COMM_WORLD,&pc));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCNONE));
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetFromOptions(pc));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetOperators(pc,mat,mat));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetUp(pc));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Create KSP context and set up data structures */
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(PETSC_COMM_WORLD,&ksp));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetType(ksp,KSPRICHARDSON));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(ksp));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetOperators(pc,mat,mat));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetPC(ksp,pc));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetUp(ksp));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Solve the problem */
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetType(ksp,&kspname));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCGetType(pc,&pcname));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolve(ksp,b,u));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(u,-1.0,ustar));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(u,NORM_2,&norm));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetIterationNumber(ksp,&its));
68c4762a1bSJed Brown   if (norm > tol) {
69*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"2 norm of error %g Number of iterations %D\n",(double)norm,its));
70c4762a1bSJed Brown   }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Free data structures */
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&ksp));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ustar));
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCDestroy(&pc));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   ierr = PetscFinalize();
81c4762a1bSJed Brown   return ierr;
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown /*TEST
85c4762a1bSJed Brown 
86c4762a1bSJed Brown    test:
87c4762a1bSJed Brown       args: -ksp_type cg -ksp_monitor_short
88c4762a1bSJed Brown 
89c4762a1bSJed Brown TEST*/
90