xref: /petsc/src/ksp/pc/tests/ex7.c (revision 856bee69f0e0908e75ff837867b1777dfb1ced96)
1 
2 static char help[] = "Tests MatILUFactorSymbolic() on matrix with missing diagonal.\n\n";
3 
4 #include <petscmat.h>
5 #include <petscpc.h>
6 
7 int main(int argc, char **args)
8 {
9   Mat         C, A;
10   PetscInt    i, j;
11   PetscScalar v;
12   PC          pc;
13   Vec         xtmp;
14 
15   PetscFunctionBeginUser;
16   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
17   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
18   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, 3, 3));
19   PetscCall(MatSetFromOptions(C));
20   PetscCall(MatSetUp(C));
21   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 3, &xtmp));
22   i = 0;
23   j = 0;
24   v = 4;
25   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
26   i = 0;
27   j = 2;
28   v = 1;
29   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
30   i = 1;
31   j = 0;
32   v = 1;
33   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
34   i = 1;
35   j = 1;
36   v = 4;
37   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
38   i = 2;
39   j = 1;
40   v = 1;
41   PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
42 
43   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
44   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
45 
46   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
47   PetscCall(PCCreate(PETSC_COMM_WORLD, &pc));
48   PetscCall(PCSetFromOptions(pc));
49   PetscCall(PCSetOperators(pc, C, C));
50   PetscCall(PCSetUp(pc));
51   PetscCall(PCFactorGetMatrix(pc, &A));
52   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
53 
54   PetscCall(PCDestroy(&pc));
55   PetscCall(VecDestroy(&xtmp));
56   PetscCall(MatDestroy(&C));
57 
58   PetscCall(PetscFinalize());
59   return 0;
60 }
61