xref: /petsc/src/mat/tests/ex22.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode    ierr;
12c4762a1bSJed Brown   PetscScalar       v;
13c4762a1bSJed Brown   PetscMPIInt       size;
14c4762a1bSJed Brown   IS                rperm,cperm,icperm;
15c4762a1bSJed Brown   const PetscInt    *rperm_ptr,*cperm_ptr,*cols;
16c4762a1bSJed Brown   const PetscScalar *vals;
17c4762a1bSJed Brown   PetscBool         TestMyorder=PETSC_FALSE;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
21*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN */
24c4762a1bSJed Brown   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
26c4762a1bSJed Brown   for (i=0; i<m; i++) {
27c4762a1bSJed Brown     for (j=0; j<n; j++) {
28c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
29c4762a1bSJed Brown       if (i>0)   {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
30c4762a1bSJed Brown       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
31c4762a1bSJed Brown       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
32c4762a1bSJed Brown       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
33c4762a1bSJed Brown       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
34c4762a1bSJed Brown     }
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   ierr = MatGetOrdering(C,MATORDERINGND,&rperm,&cperm);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = ISView(rperm,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = ISDestroy(&rperm);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   ierr = MatGetOrdering(C,MATORDERINGRCM,&rperm,&cperm);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = ISView(rperm,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = ISDestroy(&rperm);CHKERRQ(ierr);
47c4762a1bSJed Brown   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   ierr = MatGetOrdering(C,MATORDERINGQMD,&rperm,&cperm);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = ISView(rperm,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
51c4762a1bSJed Brown   ierr = ISDestroy(&rperm);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* create Cperm = rperm*C*icperm */
55c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-testmyordering",&TestMyorder,NULL);CHKERRQ(ierr);
56c4762a1bSJed Brown   if (TestMyorder) {
57c4762a1bSJed Brown     ierr = MatGetOrdering_myordering(C,MATORDERINGQMD,&rperm,&cperm);CHKERRQ(ierr);
58c4762a1bSJed Brown     printf("myordering's rperm:\n");
59c4762a1bSJed Brown     ierr = ISView(rperm,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
60c4762a1bSJed Brown     ierr = ISInvertPermutation(cperm,PETSC_DECIDE,&icperm);CHKERRQ(ierr);
61c4762a1bSJed Brown     ierr = ISGetIndices(rperm,&rperm_ptr);CHKERRQ(ierr);
62c4762a1bSJed Brown     ierr = ISGetIndices(icperm,&cperm_ptr);CHKERRQ(ierr);
63c4762a1bSJed Brown     ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&Cperm);CHKERRQ(ierr);
64c4762a1bSJed Brown     for (i=0; i<m*n; i++) {
65c4762a1bSJed Brown       ierr = MatGetRow(C,rperm_ptr[i],&ncols,&cols,&vals);CHKERRQ(ierr);
66c4762a1bSJed Brown       for (j=0; j<ncols; j++) {
67c4762a1bSJed Brown         /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */
68c4762a1bSJed Brown         ierr = MatSetValues(Cperm,1,&i,1,&cperm_ptr[cols[j]],&vals[j],INSERT_VALUES);CHKERRQ(ierr);
69c4762a1bSJed Brown       }
70c4762a1bSJed Brown     }
71c4762a1bSJed Brown     ierr = MatAssemblyBegin(Cperm,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72c4762a1bSJed Brown     ierr = MatAssemblyEnd(Cperm,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
73c4762a1bSJed Brown     ierr = ISRestoreIndices(rperm,&rperm_ptr);CHKERRQ(ierr);
74c4762a1bSJed Brown     ierr = ISRestoreIndices(icperm,&cperm_ptr);CHKERRQ(ierr);
75c4762a1bSJed Brown 
76c4762a1bSJed Brown     ierr = ISDestroy(&rperm);CHKERRQ(ierr);
77c4762a1bSJed Brown     ierr = ISDestroy(&cperm);CHKERRQ(ierr);
78c4762a1bSJed Brown     ierr = ISDestroy(&icperm);CHKERRQ(ierr);
79c4762a1bSJed Brown     ierr = MatDestroy(&Cperm);CHKERRQ(ierr);
80c4762a1bSJed Brown   }
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
83c4762a1bSJed Brown   ierr = PetscFinalize();
84c4762a1bSJed Brown   return ierr;
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown #include <petsc/private/matimpl.h>
88c4762a1bSJed Brown /* This is modified from MatGetOrdering_Natural() */
89c4762a1bSJed Brown PetscErrorCode MatGetOrdering_myordering(Mat mat,MatOrderingType type,IS *irow,IS *icol)
90c4762a1bSJed Brown {
91c4762a1bSJed Brown   PetscErrorCode ierr;
92c4762a1bSJed Brown   PetscInt       n,i,*ii;
93c4762a1bSJed Brown   PetscBool      done;
94c4762a1bSJed Brown   MPI_Comm       comm;
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   PetscFunctionBegin;
97c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,NULL,NULL,&done);CHKERRQ(ierr);
100c4762a1bSJed Brown   if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */
101c4762a1bSJed Brown     ierr = PetscMalloc1(n,&ii);CHKERRQ(ierr);
102c4762a1bSJed Brown     for (i=0; i<n; i++) ii[i] = n-i-1; /* replace your index here */
103c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_COPY_VALUES,irow);CHKERRQ(ierr);
104c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_OWN_POINTER,icol);CHKERRQ(ierr);
105c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"MatRestoreRowIJ fails!");
106c4762a1bSJed Brown   ierr = ISSetPermutation(*irow);CHKERRQ(ierr);
107c4762a1bSJed Brown   ierr = ISSetPermutation(*icol);CHKERRQ(ierr);
108c4762a1bSJed Brown   PetscFunctionReturn(0);
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111c4762a1bSJed Brown /*TEST
112c4762a1bSJed Brown 
113c4762a1bSJed Brown    test:
114c4762a1bSJed Brown 
115c4762a1bSJed Brown TEST*/
116