Actual source code: ex2.c
2: static char help[] = "Tests PC and KSP on a tridiagonal matrix. Note that most\n\
3: users should employ the KSP interface instead of using PC directly.\n\n";
5: #include petscksp.h
6: #include petscpc.h
7: #include petsc.h
8: #include <stdio.h>
12: int main(int argc,char **args)
13: {
14: Mat mat; /* matrix */
15: Vec b,ustar,u; /* vectors (RHS, exact solution, approx solution) */
16: PC pc; /* PC context */
17: KSP ksp; /* KSP context */
19: PetscInt n = 10,i,its,col[3];
20: PetscScalar value[3];
21: const PCType pcname;
22: const KSPType kspname;
23: PetscReal norm;
25: PetscInitialize(&argc,&args,(char *)0,help);
27: /* Create and initialize vectors */
28: VecCreateSeq(PETSC_COMM_SELF,n,&b);
29: VecCreateSeq(PETSC_COMM_SELF,n,&ustar);
30: VecCreateSeq(PETSC_COMM_SELF,n,&u);
31: VecSet(ustar,1.0);
32: VecSet(u,0.0);
34: /* Create and assemble matrix */
35: MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,PETSC_NULL,&mat);
36: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
37: for (i=1; i<n-1; i++) {
38: col[0] = i-1; col[1] = i; col[2] = i+1;
39: MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES);
40: }
41: i = n - 1; col[0] = n - 2; col[1] = n - 1;
42: MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);
43: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
44: MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);
45: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
46: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
48: /* Compute right-hand-side vector */
49: MatMult(mat,ustar,b);
51: /* Create PC context and set up data structures */
52: PCCreate(PETSC_COMM_WORLD,&pc);
53: PCSetType(pc,PCNONE);
54: PCSetFromOptions(pc);
55: PCSetOperators(pc,mat,mat,DIFFERENT_NONZERO_PATTERN);
56: PCSetUp(pc);
58: /* Create KSP context and set up data structures */
59: KSPCreate(PETSC_COMM_WORLD,&ksp);
60: KSPSetType(ksp,KSPRICHARDSON);
61: KSPSetFromOptions(ksp);
62: PCSetOperators(pc,mat,mat,DIFFERENT_NONZERO_PATTERN);
63: KSPSetPC(ksp,pc);
64: KSPSetUp(ksp);
66: /* Solve the problem */
67: KSPGetType(ksp,&kspname);
68: PCGetType(pc,&pcname);
69: PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname);
70: KSPSolve(ksp,b,u);
71: VecAXPY(u,-1.0,ustar);
72: VecNorm(u,NORM_2,&norm);
73: KSPGetIterationNumber(ksp,&its);
74: PetscPrintf(PETSC_COMM_SELF,"2 norm of error %A Number of iterations %D\n",norm,its);
76: /* Free data structures */
77: KSPDestroy(ksp);
78: VecDestroy(u);
79: VecDestroy(ustar);
80: VecDestroy(b);
81: MatDestroy(mat);
82: PCDestroy(pc);
84: PetscFinalize();
85: return 0;
86: }
87: