xref: /petsc/src/ksp/pc/tests/ex2.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
7*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
8*d71ae5a4SJacob Faibussowitsch {
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 
19327415f7SBarry 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));
309371c9d4SSatish Balay   value[0] = -1.0;
319371c9d4SSatish Balay   value[1] = 2.0;
329371c9d4SSatish Balay   value[2] = -1.0;
33c4762a1bSJed Brown   for (i = 1; i < n - 1; i++) {
349371c9d4SSatish Balay     col[0] = i - 1;
359371c9d4SSatish Balay     col[1] = i;
369371c9d4SSatish Balay     col[2] = i + 1;
379566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat, 1, &i, 3, col, value, INSERT_VALUES));
38c4762a1bSJed Brown   }
399371c9d4SSatish Balay   i      = n - 1;
409371c9d4SSatish Balay   col[0] = n - 2;
419371c9d4SSatish Balay   col[1] = n - 1;
429566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES));
439371c9d4SSatish Balay   i        = 0;
449371c9d4SSatish Balay   col[0]   = 0;
459371c9d4SSatish Balay   col[1]   = 1;
469371c9d4SSatish Balay   value[0] = 2.0;
479371c9d4SSatish Balay   value[1] = -1.0;
489566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES));
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Compute right-hand-side vector */
539566063dSJacob Faibussowitsch   PetscCall(MatMult(mat, ustar, b));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Create PC context and set up data structures */
569566063dSJacob Faibussowitsch   PetscCall(PCCreate(PETSC_COMM_WORLD, &pc));
579566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
589566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
599566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc, mat, mat));
609566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* Create KSP context and set up data structures */
639566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
649566063dSJacob Faibussowitsch   PetscCall(KSPSetType(ksp, KSPRICHARDSON));
659566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
669566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc, mat, mat));
679566063dSJacob Faibussowitsch   PetscCall(KSPSetPC(ksp, pc));
689566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(ksp));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* Solve the problem */
719566063dSJacob Faibussowitsch   PetscCall(KSPGetType(ksp, &kspname));
729566063dSJacob Faibussowitsch   PetscCall(PCGetType(pc, &pcname));
739566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Running %s with %s preconditioning\n", kspname, pcname));
749566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, b, u));
759566063dSJacob Faibussowitsch   PetscCall(VecAXPY(u, -1.0, ustar));
769566063dSJacob Faibussowitsch   PetscCall(VecNorm(u, NORM_2, &norm));
779566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(ksp, &its));
7848a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "2 norm of error %g Number of iterations %" PetscInt_FMT "\n", (double)norm, its));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Free data structures */
819566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ustar));
849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
869566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
89b122ec5aSJacob Faibussowitsch   return 0;
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
92c4762a1bSJed Brown /*TEST
93c4762a1bSJed Brown 
94c4762a1bSJed Brown    test:
95c4762a1bSJed Brown       args: -ksp_type cg -ksp_monitor_short
96c4762a1bSJed Brown 
97c4762a1bSJed Brown TEST*/
98