xref: /petsc/src/ksp/pc/tests/ex7.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatILUFactorSymbolic() on matrix with missing diagonal.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown #include <petscpc.h>
6c4762a1bSJed Brown 
7*9371c9d4SSatish Balay int main(int argc, char **args) {
8c4762a1bSJed Brown   Mat         C, A;
9c4762a1bSJed Brown   PetscInt    i, j;
10c4762a1bSJed Brown   PetscScalar v;
11c4762a1bSJed Brown   PC          pc;
12c4762a1bSJed Brown   Vec         xtmp;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
169566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
179566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, 3, 3));
189566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
199566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
209566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 3, &xtmp));
21*9371c9d4SSatish Balay   i = 0;
22*9371c9d4SSatish Balay   j = 0;
23*9371c9d4SSatish Balay   v = 4;
249566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
25*9371c9d4SSatish Balay   i = 0;
26*9371c9d4SSatish Balay   j = 2;
27*9371c9d4SSatish Balay   v = 1;
289566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
29*9371c9d4SSatish Balay   i = 1;
30*9371c9d4SSatish Balay   j = 0;
31*9371c9d4SSatish Balay   v = 1;
329566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
33*9371c9d4SSatish Balay   i = 1;
34*9371c9d4SSatish Balay   j = 1;
35*9371c9d4SSatish Balay   v = 4;
369566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
37*9371c9d4SSatish Balay   i = 2;
38*9371c9d4SSatish Balay   j = 1;
39*9371c9d4SSatish Balay   v = 1;
409566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
469566063dSJacob Faibussowitsch   PetscCall(PCCreate(PETSC_COMM_WORLD, &pc));
479566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
489566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc, C, C));
499566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
509566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatrix(pc, &A));
519566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc));
549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xtmp));
559566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
58b122ec5aSJacob Faibussowitsch   return 0;
59c4762a1bSJed Brown }
60