xref: /petsc/src/ksp/pc/tests/ex7.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,3,3));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD,3,&xtmp));
21c4762a1bSJed Brown   i    = 0; j = 0; v = 4;
225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
23c4762a1bSJed Brown   i    = 0; j = 2; v = 1;
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
25c4762a1bSJed Brown   i    = 1; j = 0; v = 1;
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
27c4762a1bSJed Brown   i    = 1; j = 1; v = 4;
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
29c4762a1bSJed Brown   i    = 2; j = 1; v = 1;
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES));
31c4762a1bSJed Brown 
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
34c4762a1bSJed Brown 
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(PCCreate(PETSC_COMM_WORLD,&pc));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetFromOptions(pc));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetOperators(pc,C,C));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetUp(pc));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PCFactorGetMatrix(pc,&A));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
42c4762a1bSJed Brown 
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PCDestroy(&pc));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&xtmp));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
46c4762a1bSJed Brown 
47*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
48*b122ec5aSJacob Faibussowitsch   return 0;
49c4762a1bSJed Brown }
50