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