xref: /petsc/src/mat/tests/ex6.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests reordering a matrix.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            C;
9c4762a1bSJed Brown   PetscInt       i,j,m = 5,n = 5,Ii,J;
10c4762a1bSJed Brown   PetscScalar    v;
11c4762a1bSJed Brown   IS             perm,iperm;
12c4762a1bSJed Brown 
13*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
20c4762a1bSJed Brown   for (i=0; i<m; i++) {
21c4762a1bSJed Brown     for (j=0; j<n; j++) {
22c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
235f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
245f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
255f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
265f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
275f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
28c4762a1bSJed Brown     }
29c4762a1bSJed Brown   }
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
32c4762a1bSJed Brown 
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(C,MATORDERINGND,&perm,&iperm));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(perm,PETSC_VIEWER_STDOUT_SELF));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(iperm,PETSC_VIEWER_STDOUT_SELF));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_SELF));
37c4762a1bSJed Brown 
385f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&perm));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iperm));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
41*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
42*b122ec5aSJacob Faibussowitsch   return 0;
43c4762a1bSJed Brown }
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /*TEST
46c4762a1bSJed Brown 
47c4762a1bSJed Brown    test:
48c4762a1bSJed Brown 
49c4762a1bSJed Brown TEST*/
50