Actual source code: ex2.c
1: /*$Id: ex2.c,v 1.53 2001/08/07 21:30:37 bsmith Exp $*/
3: static char help[] = "Tests PC and KSP on a tridiagonal matrix. Note that mostn
4: users should employ the SLES interface instead of using PC directly.nn";
6: #include petscksp.h
7: #include petscpc.h
8: #include petsc.h
9: #include <stdio.h>
11: int main(int argc,char **args)
12: {
13: Mat mat; /* matrix */
14: Vec b,ustar,u; /* vectors (RHS, exact solution, approx solution) */
15: PC pc; /* PC context */
16: KSP ksp; /* KSP context */
17: int ierr,n = 10,i,its,col[3];
18: PetscScalar value[3],mone = -1.0,one = 1.0,zero = 0.0;
19: PCType pcname;
20: KSPType kspname;
21: PetscReal norm;
23: PetscInitialize(&argc,&args,(char *)0,help);
25: /* Create and initialize vectors */
26: VecCreateSeq(PETSC_COMM_SELF,n,&b);
27: VecCreateSeq(PETSC_COMM_SELF,n,&ustar);
28: VecCreateSeq(PETSC_COMM_SELF,n,&u);
29: VecSet(&one,ustar);
30: VecSet(&zero,u);
32: /* Create and assemble matrix */
33: MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,PETSC_NULL,&mat);
34: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
35: for (i=1; i<n-1; i++) {
36: col[0] = i-1; col[1] = i; col[2] = i+1;
37: MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES);
38: }
39: i = n - 1; col[0] = n - 2; col[1] = n - 1;
40: MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);
41: i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
42: MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);
43: MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
44: MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
46: /* Compute right-hand-side vector */
47: MatMult(mat,ustar,b);
49: /* Create PC context and set up data structures */
50: PCCreate(PETSC_COMM_WORLD,&pc);
51: PCSetType(pc,PCNONE);
52: PCSetFromOptions(pc);
53: PCSetOperators(pc,mat,mat,DIFFERENT_NONZERO_PATTERN);
54: PCSetVector(pc,u);
55: PCSetUp(pc);
57: /* Create KSP context and set up data structures */
58: KSPCreate(PETSC_COMM_WORLD,&ksp);
59: KSPSetType(ksp,KSPRICHARDSON);
60: KSPSetFromOptions(ksp);
61: KSPSetSolution(ksp,u);
62: KSPSetRhs(ksp,b);
63: PCSetOperators(pc,mat,mat,DIFFERENT_NONZERO_PATTERN);
64: KSPSetPC(ksp,pc);
65: KSPSetUp(ksp);
67: /* Solve the problem */
68: KSPGetType(ksp,&kspname);
69: PCGetType(pc,&pcname);
70: PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioningn",kspname,pcname);
71: KSPSolve(ksp,&its);
72: VecAXPY(&mone,ustar,u);
73: VecNorm(u,NORM_2,&norm);
74: PetscPrintf(PETSC_COMM_SELF,"2 norm of error %A Number of iterations %dn",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: