xref: /petsc/src/ksp/pc/tests/ex7.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
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   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
16   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
17   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,3,3));
18   PetscCall(MatSetFromOptions(C));
19   PetscCall(MatSetUp(C));
20   PetscCall(VecCreateSeq(PETSC_COMM_WORLD,3,&xtmp));
21   i    = 0; j = 0; v = 4;
22   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
23   i    = 0; j = 2; v = 1;
24   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
25   i    = 1; j = 0; v = 1;
26   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
27   i    = 1; j = 1; v = 4;
28   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
29   i    = 2; j = 1; v = 1;
30   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
31 
32   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
33   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
34 
35   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
36   PetscCall(PCCreate(PETSC_COMM_WORLD,&pc));
37   PetscCall(PCSetFromOptions(pc));
38   PetscCall(PCSetOperators(pc,C,C));
39   PetscCall(PCSetUp(pc));
40   PetscCall(PCFactorGetMatrix(pc,&A));
41   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
42 
43   PetscCall(PCDestroy(&pc));
44   PetscCall(VecDestroy(&xtmp));
45   PetscCall(MatDestroy(&C));
46 
47   PetscCall(PetscFinalize());
48   return 0;
49 }
50