xref: /petsc/src/ksp/pc/tests/ex7.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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; j = 0; v = 4;
23   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
24   i    = 0; j = 2; v = 1;
25   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
26   i    = 1; j = 0; v = 1;
27   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
28   i    = 1; j = 1; v = 4;
29   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
30   i    = 2; j = 1; v = 1;
31   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
32 
33   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
34   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
35 
36   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
37   PetscCall(PCCreate(PETSC_COMM_WORLD,&pc));
38   PetscCall(PCSetFromOptions(pc));
39   PetscCall(PCSetOperators(pc,C,C));
40   PetscCall(PCSetUp(pc));
41   PetscCall(PCFactorGetMatrix(pc,&A));
42   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
43 
44   PetscCall(PCDestroy(&pc));
45   PetscCall(VecDestroy(&xtmp));
46   PetscCall(MatDestroy(&C));
47 
48   PetscCall(PetscFinalize());
49   return 0;
50 }
51