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