xref: /petsc/src/ksp/pc/tests/ex7.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
8*d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   Mat         C, A;
10c4762a1bSJed Brown   PetscInt    i, j;
11c4762a1bSJed Brown   PetscScalar v;
12c4762a1bSJed Brown   PC          pc;
13c4762a1bSJed Brown   Vec         xtmp;
14c4762a1bSJed Brown 
15327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, 3, 3));
199566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
209566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
219566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 3, &xtmp));
229371c9d4SSatish Balay   i = 0;
239371c9d4SSatish Balay   j = 0;
249371c9d4SSatish Balay   v = 4;
259566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
269371c9d4SSatish Balay   i = 0;
279371c9d4SSatish Balay   j = 2;
289371c9d4SSatish Balay   v = 1;
299566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
309371c9d4SSatish Balay   i = 1;
319371c9d4SSatish Balay   j = 0;
329371c9d4SSatish Balay   v = 1;
339566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
349371c9d4SSatish Balay   i = 1;
359371c9d4SSatish Balay   j = 1;
369371c9d4SSatish Balay   v = 4;
379566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
389371c9d4SSatish Balay   i = 2;
399371c9d4SSatish Balay   j = 1;
409371c9d4SSatish Balay   v = 1;
419566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
42c4762a1bSJed Brown 
439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
479566063dSJacob Faibussowitsch   PetscCall(PCCreate(PETSC_COMM_WORLD, &pc));
489566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
499566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc, C, C));
509566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
519566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatrix(pc, &A));
529566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc));
559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xtmp));
569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
59b122ec5aSJacob Faibussowitsch   return 0;
60c4762a1bSJed Brown }
61