xref: /petsc/src/mat/tests/ex22.c (revision be096a4675ce3b20dd7f9c195e6046e69d9955cf)
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 
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
20*be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN */
239566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C));
249566063dSJacob Faibussowitsch   PetscCall(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;
289566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
299566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
309566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
319566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
329566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
33c4762a1bSJed Brown     }
34c4762a1bSJed Brown   }
359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
37c4762a1bSJed Brown 
389566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C,MATORDERINGND,&rperm,&cperm));
399566063dSJacob Faibussowitsch   PetscCall(ISView(rperm,PETSC_VIEWER_STDOUT_SELF));
409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&rperm));
419566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
42c4762a1bSJed Brown 
439566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C,MATORDERINGRCM,&rperm,&cperm));
449566063dSJacob Faibussowitsch   PetscCall(ISView(rperm,PETSC_VIEWER_STDOUT_SELF));
459566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&rperm));
469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C,MATORDERINGQMD,&rperm,&cperm));
499566063dSJacob Faibussowitsch   PetscCall(ISView(rperm,PETSC_VIEWER_STDOUT_SELF));
509566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&rperm));
519566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* create Cperm = rperm*C*icperm */
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-testmyordering",&TestMyorder,NULL));
55c4762a1bSJed Brown   if (TestMyorder) {
569566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering_myordering(C,MATORDERINGQMD,&rperm,&cperm));
57c4762a1bSJed Brown     printf("myordering's rperm:\n");
589566063dSJacob Faibussowitsch     PetscCall(ISView(rperm,PETSC_VIEWER_STDOUT_SELF));
599566063dSJacob Faibussowitsch     PetscCall(ISInvertPermutation(cperm,PETSC_DECIDE,&icperm));
609566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(rperm,&rperm_ptr));
619566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icperm,&cperm_ptr));
629566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&Cperm));
63c4762a1bSJed Brown     for (i=0; i<m*n; i++) {
649566063dSJacob Faibussowitsch       PetscCall(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]); */
679566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Cperm,1,&i,1,&cperm_ptr[cols[j]],&vals[j],INSERT_VALUES));
68c4762a1bSJed Brown       }
69c4762a1bSJed Brown     }
709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Cperm,MAT_FINAL_ASSEMBLY));
719566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Cperm,MAT_FINAL_ASSEMBLY));
729566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(rperm,&rperm_ptr));
739566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icperm,&cperm_ptr));
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&rperm));
769566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cperm));
779566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&icperm));
789566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Cperm));
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
829566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
83b122ec5aSJacob 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;
959566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)mat,&comm));
969566063dSJacob Faibussowitsch   PetscCall(MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done));
979566063dSJacob Faibussowitsch   PetscCall(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 */
999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&ii));
100c4762a1bSJed Brown     for (i=0; i<n; i++) ii[i] = n-i-1; /* replace your index here */
1019566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_COPY_VALUES,irow));
1029566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_OWN_POINTER,icol));
103c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"MatRestoreRowIJ fails!");
1049566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*irow));
1059566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(*icol));
106c4762a1bSJed Brown   PetscFunctionReturn(0);
107c4762a1bSJed Brown }
108c4762a1bSJed Brown 
109c4762a1bSJed Brown /*TEST
110c4762a1bSJed Brown 
111c4762a1bSJed Brown    test:
112c4762a1bSJed Brown 
113c4762a1bSJed Brown TEST*/
114