xref: /petsc/src/ksp/pc/tests/ex7.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Tests MatILUFactorSymbolic() on matrix with missing diagonal.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
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         C, A;
9c4762a1bSJed Brown   PetscInt    i, j;
10c4762a1bSJed Brown   PetscScalar v;
11c4762a1bSJed Brown   PC          pc;
12c4762a1bSJed Brown   Vec         xtmp;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
15*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, 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));
219371c9d4SSatish Balay   i = 0;
229371c9d4SSatish Balay   j = 0;
239371c9d4SSatish Balay   v = 4;
249566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
259371c9d4SSatish Balay   i = 0;
269371c9d4SSatish Balay   j = 2;
279371c9d4SSatish Balay   v = 1;
289566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
299371c9d4SSatish Balay   i = 1;
309371c9d4SSatish Balay   j = 0;
319371c9d4SSatish Balay   v = 1;
329566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
339371c9d4SSatish Balay   i = 1;
349371c9d4SSatish Balay   j = 1;
359371c9d4SSatish Balay   v = 4;
369566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
379371c9d4SSatish Balay   i = 2;
389371c9d4SSatish Balay   j = 1;
399371c9d4SSatish 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