xref: /petsc/src/ksp/pc/tests/ex2.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
21c4762a1bSJed Brown   /* Create and initialize vectors */
229566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&b));
239566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&ustar));
249566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&u));
259566063dSJacob Faibussowitsch   PetscCall(VecSet(ustar,1.0));
269566063dSJacob Faibussowitsch   PetscCall(VecSet(u,0.0));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Create and assemble matrix */
299566063dSJacob Faibussowitsch   PetscCall(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;
339566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES));
34c4762a1bSJed Brown   }
35c4762a1bSJed Brown   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
369566063dSJacob Faibussowitsch   PetscCall(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;
389566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES));
399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* Compute right-hand-side vector */
439566063dSJacob Faibussowitsch   PetscCall(MatMult(mat,ustar,b));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Create PC context and set up data structures */
469566063dSJacob Faibussowitsch   PetscCall(PCCreate(PETSC_COMM_WORLD,&pc));
479566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCNONE));
489566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
499566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc,mat,mat));
509566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Create KSP context and set up data structures */
539566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp));
549566063dSJacob Faibussowitsch   PetscCall(KSPSetType(ksp,KSPRICHARDSON));
559566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
569566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc,mat,mat));
579566063dSJacob Faibussowitsch   PetscCall(KSPSetPC(ksp,pc));
589566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(ksp));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Solve the problem */
619566063dSJacob Faibussowitsch   PetscCall(KSPGetType(ksp,&kspname));
629566063dSJacob Faibussowitsch   PetscCall(PCGetType(pc,&pcname));
639566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname));
649566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp,b,u));
659566063dSJacob Faibussowitsch   PetscCall(VecAXPY(u,-1.0,ustar));
669566063dSJacob Faibussowitsch   PetscCall(VecNorm(u,NORM_2,&norm));
679566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(ksp,&its));
68c4762a1bSJed Brown   if (norm > tol) {
6963a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"2 norm of error %g Number of iterations %" PetscInt_FMT "\n",(double)norm,its));
70c4762a1bSJed Brown   }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Free data structures */
739566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ustar));
769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
789566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc));
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
81b122ec5aSJacob Faibussowitsch   return 0;
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