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