xref: /petsc/src/mat/tests/ex22.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests matrix ordering routines.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown extern PetscErrorCode MatGetOrdering_myordering(Mat,MatOrderingType,IS*,IS*);
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat               C,Cperm;
10c4762a1bSJed Brown   PetscInt          i,j,m = 5,n = 5,Ii,J,ncols;
11c4762a1bSJed Brown   PetscScalar       v;
12c4762a1bSJed Brown   PetscMPIInt       size;
13c4762a1bSJed Brown   IS                rperm,cperm,icperm;
14c4762a1bSJed Brown   const PetscInt    *rperm_ptr,*cperm_ptr,*cols;
15c4762a1bSJed Brown   const PetscScalar *vals;
16c4762a1bSJed Brown   PetscBool         TestMyorder=PETSC_FALSE;
17c4762a1bSJed Brown 
18*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
195f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
202c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN */
235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
25c4762a1bSJed Brown   for (i=0; i<m; i++) {
26c4762a1bSJed Brown     for (j=0; j<n; j++) {
27c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
285f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
295f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
305f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
315f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
325f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
33c4762a1bSJed Brown     }
34c4762a1bSJed Brown   }
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
37c4762a1bSJed Brown 
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(C,MATORDERINGND,&rperm,&cperm));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&rperm));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&cperm));
42c4762a1bSJed Brown 
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(C,MATORDERINGRCM,&rperm,&cperm));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&rperm));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&cperm));
47c4762a1bSJed Brown 
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(C,MATORDERINGQMD,&rperm,&cperm));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&rperm));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&cperm));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* create Cperm = rperm*C*icperm */
545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-testmyordering",&TestMyorder,NULL));
55c4762a1bSJed Brown   if (TestMyorder) {
565f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetOrdering_myordering(C,MATORDERINGQMD,&rperm,&cperm));
57c4762a1bSJed Brown     printf("myordering's rperm:\n");
585f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF));
595f80ce2aSJacob Faibussowitsch     CHKERRQ(ISInvertPermutation(cperm,PETSC_DECIDE,&icperm));
605f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(rperm,&rperm_ptr));
615f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(icperm,&cperm_ptr));
625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&Cperm));
63c4762a1bSJed Brown     for (i=0; i<m*n; i++) {
645f80ce2aSJacob Faibussowitsch       CHKERRQ(MatGetRow(C,rperm_ptr[i],&ncols,&cols,&vals));
65c4762a1bSJed Brown       for (j=0; j<ncols; j++) {
66c4762a1bSJed Brown         /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */
675f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(Cperm,1,&i,1,&cperm_ptr[cols[j]],&vals[j],INSERT_VALUES));
68c4762a1bSJed Brown       }
69c4762a1bSJed Brown     }
705f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(Cperm,MAT_FINAL_ASSEMBLY));
715f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(Cperm,MAT_FINAL_ASSEMBLY));
725f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(rperm,&rperm_ptr));
735f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(icperm,&cperm_ptr));
74c4762a1bSJed Brown 
755f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&rperm));
765f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&cperm));
775f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&icperm));
785f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Cperm));
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown 
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
82*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
83*b122ec5aSJacob Faibussowitsch   return 0;
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown #include <petsc/private/matimpl.h>
87c4762a1bSJed Brown /* This is modified from MatGetOrdering_Natural() */
88c4762a1bSJed Brown PetscErrorCode MatGetOrdering_myordering(Mat mat,MatOrderingType type,IS *irow,IS *icol)
89c4762a1bSJed Brown {
90c4762a1bSJed Brown   PetscInt       n,i,*ii;
91c4762a1bSJed Brown   PetscBool      done;
92c4762a1bSJed Brown   MPI_Comm       comm;
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   PetscFunctionBegin;
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)mat,&comm));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,NULL,NULL,&done));
98c4762a1bSJed Brown   if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n,&ii));
100c4762a1bSJed Brown     for (i=0; i<n; i++) ii[i] = n-i-1; /* replace your index here */
1015f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_COPY_VALUES,irow));
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_OWN_POINTER,icol));
103c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"MatRestoreRowIJ fails!");
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(ISSetPermutation(*irow));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(ISSetPermutation(*icol));
106c4762a1bSJed Brown   PetscFunctionReturn(0);
107c4762a1bSJed Brown }
108c4762a1bSJed Brown 
109c4762a1bSJed Brown /*TEST
110c4762a1bSJed Brown 
111c4762a1bSJed Brown    test:
112c4762a1bSJed Brown 
113c4762a1bSJed Brown TEST*/
114