xref: /petsc/src/ksp/pc/tests/ex7.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            C,A;
10c4762a1bSJed Brown   PetscInt       i,j;
11c4762a1bSJed Brown   PetscScalar    v;
12c4762a1bSJed Brown   PC             pc;
13c4762a1bSJed Brown   Vec            xtmp;
14c4762a1bSJed Brown 
15*327415f7SBarry 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));
22c4762a1bSJed Brown   i    = 0; j = 0; v = 4;
239566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
24c4762a1bSJed Brown   i    = 0; j = 2; v = 1;
259566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
26c4762a1bSJed Brown   i    = 1; j = 0; v = 1;
279566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
28c4762a1bSJed Brown   i    = 1; j = 1; v = 4;
299566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
30c4762a1bSJed Brown   i    = 2; j = 1; v = 1;
319566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
35c4762a1bSJed Brown 
369566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
379566063dSJacob Faibussowitsch   PetscCall(PCCreate(PETSC_COMM_WORLD,&pc));
389566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(pc));
399566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(pc,C,C));
409566063dSJacob Faibussowitsch   PetscCall(PCSetUp(pc));
419566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatrix(pc,&A));
429566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
43c4762a1bSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc));
459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xtmp));
469566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
49b122ec5aSJacob Faibussowitsch   return 0;
50c4762a1bSJed Brown }
51