xref: /petsc/src/ksp/pc/tests/ex4.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Demonstrates the use of fast Richardson for SOR. And tests\n\
2c4762a1bSJed Brown the MatSOR() routines.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscpc.h>
5c4762a1bSJed Brown 
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         mat;
9c4762a1bSJed Brown   Vec         b, u;
10c4762a1bSJed Brown   PC          pc;
11c4762a1bSJed Brown   PetscInt    n = 5, i, col[3];
12c4762a1bSJed Brown   PetscScalar value[3];
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
15*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
16c4762a1bSJed Brown   /* Create vectors */
179566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &b));
189566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &u));
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /* Create and assemble matrix */
219566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, n, n, NULL, &mat));
229371c9d4SSatish Balay   value[0] = -1.0;
239371c9d4SSatish Balay   value[1] = 2.0;
249371c9d4SSatish Balay   value[2] = -1.0;
25c4762a1bSJed Brown   for (i = 1; i < n - 1; i++) {
269371c9d4SSatish Balay     col[0] = i - 1;
279371c9d4SSatish Balay     col[1] = i;
289371c9d4SSatish Balay     col[2] = i + 1;
299566063dSJacob Faibussowitsch     PetscCall(MatSetValues(mat, 1, &i, 3, col, value, INSERT_VALUES));
30c4762a1bSJed Brown   }
319371c9d4SSatish Balay   i      = n - 1;
329371c9d4SSatish Balay   col[0] = n - 2;
339371c9d4SSatish Balay   col[1] = n - 1;
349566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES));
359371c9d4SSatish Balay   i        = 0;
369371c9d4SSatish Balay   col[0]   = 0;
379371c9d4SSatish Balay   col[1]   = 1;
389371c9d4SSatish Balay   value[0] = 2.0;
399371c9d4SSatish Balay   value[1] = -1.0;
409566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 1, &i, 2, col, value, INSERT_VALUES));
419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Create PC context and set up data structures */
459566063dSJacob Faibussowitsch   PetscCall(PCCreate(PETSC_COMM_WORLD, &pc));
469566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCSOR));
479566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
489566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc, mat, mat));
499566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   value[0] = 1.0;
52c4762a1bSJed Brown   for (i = 0; i < n; i++) {
539566063dSJacob Faibussowitsch     PetscCall(VecSet(u, 0.0));
549566063dSJacob Faibussowitsch     PetscCall(VecSetValues(u, 1, &i, value, INSERT_VALUES));
559566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(u));
569566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(u));
579566063dSJacob Faibussowitsch     PetscCall(PCApply(pc, u, b));
589566063dSJacob Faibussowitsch     PetscCall(VecView(b, PETSC_VIEWER_STDOUT_SELF));
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* Free data structures */
629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
639566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc));
649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
669566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
67b122ec5aSJacob Faibussowitsch   return 0;
68c4762a1bSJed Brown }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown /*TEST
71c4762a1bSJed Brown 
72c4762a1bSJed Brown    test:
73c4762a1bSJed Brown 
74c4762a1bSJed Brown TEST*/
75