xref: /petsc/src/ksp/pc/tests/ex3.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Demonstrates the use of fast Richardson for SOR. And\n\
3c4762a1bSJed Brown also tests the MatSOR() routines.  Input parameters are:\n\
4c4762a1bSJed Brown  -n <n> : problem dimension\n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #include <petscksp.h>
7c4762a1bSJed Brown #include <petscpc.h>
8c4762a1bSJed Brown 
9*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
10*d71ae5a4SJacob Faibussowitsch {
11c4762a1bSJed Brown   Mat         mat;         /* matrix */
12c4762a1bSJed Brown   Vec         b, ustar, u; /* vectors (RHS, exact solution, approx solution) */
13c4762a1bSJed Brown   PC          pc;          /* PC context */
14c4762a1bSJed Brown   KSP         ksp;         /* KSP context */
15c4762a1bSJed Brown   PetscInt    n = 10, i, its, col[3];
16c4762a1bSJed Brown   PetscScalar value[3];
17c4762a1bSJed Brown   KSPType     kspname;
18c4762a1bSJed Brown   PCType      pcname;
19c4762a1bSJed Brown 
20327415f7SBarry Smith   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* Create and initialize vectors */
259566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &b));
269566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &ustar));
279566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &u));
289566063dSJacob Faibussowitsch   PetscCall(VecSet(ustar, 1.0));
299566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* Create and assemble matrix */
329566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &mat));
339566063dSJacob Faibussowitsch   PetscCall(MatSetType(mat, MATSEQAIJ));
349566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(mat, n, n, n, n));
359566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(mat));
369566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(mat, 3, NULL));
379566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(mat, 1, 3, NULL));
389566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(mat, 1, 3, NULL));
399566063dSJacob Faibussowitsch   PetscCall(MatSeqSELLSetPreallocation(mat, 3, NULL));
409371c9d4SSatish Balay   value[0] = -1.0;
419371c9d4SSatish Balay   value[1] = 2.0;
429371c9d4SSatish Balay   value[2] = -1.0;
43c4762a1bSJed Brown   for (i = 1; i < n - 1; i++) {
449371c9d4SSatish Balay     col[0] = i - 1;
459371c9d4SSatish Balay     col[1] = i;
469371c9d4SSatish Balay     col[2] = i + 1;
479566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat, 1, &i, 3, col, value, INSERT_VALUES));
48c4762a1bSJed Brown   }
499371c9d4SSatish Balay   i      = n - 1;
509371c9d4SSatish Balay   col[0] = n - 2;
519371c9d4SSatish Balay   col[1] = n - 1;
529566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES));
539371c9d4SSatish Balay   i        = 0;
549371c9d4SSatish Balay   col[0]   = 0;
559371c9d4SSatish Balay   col[1]   = 1;
569371c9d4SSatish Balay   value[0] = 2.0;
579371c9d4SSatish Balay   value[1] = -1.0;
589566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES));
599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* Compute right-hand-side vector */
639566063dSJacob Faibussowitsch   PetscCall(MatMult(mat, ustar, b));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Create PC context and set up data structures */
669566063dSJacob Faibussowitsch   PetscCall(PCCreate(PETSC_COMM_WORLD, &pc));
679566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCNONE));
689566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
699566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc, mat, mat));
709566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Create KSP context and set up data structures */
739566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
749566063dSJacob Faibussowitsch   PetscCall(KSPSetType(ksp, KSPRICHARDSON));
759566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
769566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc, mat, mat));
779566063dSJacob Faibussowitsch   PetscCall(KSPSetPC(ksp, pc));
789566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(ksp));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Solve the problem */
819566063dSJacob Faibussowitsch   PetscCall(KSPGetType(ksp, &kspname));
829566063dSJacob Faibussowitsch   PetscCall(PCGetType(pc, &pcname));
839566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Running %s with %s preconditioning\n", kspname, pcname));
849566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, b, u));
859566063dSJacob Faibussowitsch   PetscCall(KSPGetIterationNumber(ksp, &its));
8663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of iterations %" PetscInt_FMT "\n", its));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* Free data structures */
899566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ustar));
929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
949566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc));
959566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
96b122ec5aSJacob Faibussowitsch   return 0;
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown /*TEST
100c4762a1bSJed Brown 
101c4762a1bSJed Brown    testset:
102c4762a1bSJed Brown      args: -ksp_type gmres -ksp_monitor_short -pc_type sor -pc_sor_symmetric
103c4762a1bSJed Brown      output_file: output/ex3_1.out
104c4762a1bSJed Brown      test:
105c4762a1bSJed Brown        suffix: sor_aij
106c4762a1bSJed Brown      test:
107c4762a1bSJed Brown        suffix: sor_seqbaij
108c4762a1bSJed Brown        args: -mat_type seqbaij
109c4762a1bSJed Brown      test:
110c4762a1bSJed Brown        suffix: sor_seqsbaij
111c4762a1bSJed Brown        args: -mat_type seqbaij
112c4762a1bSJed Brown      test:
113c4762a1bSJed Brown        suffix: sor_seqsell
114c4762a1bSJed Brown        args: -mat_type seqsell
115c4762a1bSJed Brown 
116c4762a1bSJed Brown TEST*/
117