xref: /petsc/src/mat/tests/ex167.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Extract submatrices using unsorted indices. For SEQSBAIJ either sort both rows and columns, or sort none.\n\n";
3*c4762a1bSJed Brown /*
4*c4762a1bSJed Brown    Take a 4x4 grid and form a 5-point stencil graph Laplacian over it.
5*c4762a1bSJed Brown    Partition the grid into two subdomains by splitting into two in the j-direction (slowest varying).
6*c4762a1bSJed Brown    Impose an overlap of 1 and order the subdomains with the j-direction varying fastest.
7*c4762a1bSJed Brown    Extract the subdomain submatrices, one per rank.
8*c4762a1bSJed Brown */
9*c4762a1bSJed Brown /* Results:
10*c4762a1bSJed Brown     Sequential:
11*c4762a1bSJed Brown     - seqaij:   will error out, if rows or columns are unsorted
12*c4762a1bSJed Brown     - seqbaij:  will extract submatrices correctly even for unsorted row or column indices
13*c4762a1bSJed Brown     - seqsbaij: will extract submatrices correctly even for unsorted row and column indices (both must be sorted or not);
14*c4762a1bSJed Brown                 CANNOT automatically report inversions, because MatGetRow is not available.
15*c4762a1bSJed Brown     MPI:
16*c4762a1bSJed Brown     - mpiaij:   will error out, if columns are unsorted
17*c4762a1bSJed Brown     - mpibaij:  will error out, if columns are unsorted.
18*c4762a1bSJed Brown     - mpisbaij: will error out, if columns are unsorted; even with unsorted rows will produce correct submatrices;
19*c4762a1bSJed Brown                 CANNOT automatically report inversions, because MatGetRow is not available.
20*c4762a1bSJed Brown */
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown #include <petscmat.h>
23*c4762a1bSJed Brown #include <petscis.h>
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown int main(int argc,char **args)
26*c4762a1bSJed Brown {
27*c4762a1bSJed Brown   Mat            A, *S;
28*c4762a1bSJed Brown   IS             rowis[2], colis[2];
29*c4762a1bSJed Brown   PetscInt       n,N,i,j,k,l,nsub,Jlow[2] = {0,1}, *jlow, Jhigh[2] = {3,4}, *jhigh, row, col, *subindices, ncols;
30*c4762a1bSJed Brown   const PetscInt *cols;
31*c4762a1bSJed Brown   PetscScalar    v;
32*c4762a1bSJed Brown   PetscMPIInt    rank, size, p, inversions, total_inversions;
33*c4762a1bSJed Brown   PetscBool      sort_rows, sort_cols, show_inversions;
34*c4762a1bSJed Brown   PetscErrorCode ierr;
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
37*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
39*c4762a1bSJed Brown   if (size>2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG, "A uniprocessor or two-processor example only.\n");
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
42*c4762a1bSJed Brown   if (size > 1) {
43*c4762a1bSJed Brown     n = 8; N = 16;
44*c4762a1bSJed Brown   } else {
45*c4762a1bSJed Brown     n = 16; N = 16;
46*c4762a1bSJed Brown   }
47*c4762a1bSJed Brown   ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   /* Don't care if the entries are set multiple times by different procs. */
52*c4762a1bSJed Brown   for (i=0; i<4; ++i) {
53*c4762a1bSJed Brown     for (j = 0; j<4; ++j) {
54*c4762a1bSJed Brown       row = j*4+i;
55*c4762a1bSJed Brown       v   = -1.0;
56*c4762a1bSJed Brown       if (i>0) {
57*c4762a1bSJed Brown         col =  row-1; ierr = MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr);
58*c4762a1bSJed Brown       }
59*c4762a1bSJed Brown       if (i<3) {
60*c4762a1bSJed Brown         col = row+1; ierr = MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr);
61*c4762a1bSJed Brown       }
62*c4762a1bSJed Brown       if (j>0) {
63*c4762a1bSJed Brown         col = row-4; ierr = MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr);
64*c4762a1bSJed Brown       }
65*c4762a1bSJed Brown       if (j<3) {
66*c4762a1bSJed Brown         col = row+4; ierr = MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr);
67*c4762a1bSJed Brown       }
68*c4762a1bSJed Brown       v    = 4.0;
69*c4762a1bSJed Brown       ierr = MatSetValues(A,1,&row,1,&row,&v,INSERT_VALUES);CHKERRQ(ierr);
70*c4762a1bSJed Brown     }
71*c4762a1bSJed Brown   }
72*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "Original matrix\n");CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   if (size > 1) {
78*c4762a1bSJed Brown     nsub = 1; /* one subdomain per rank */
79*c4762a1bSJed Brown   }
80*c4762a1bSJed Brown   else {
81*c4762a1bSJed Brown     nsub = 2; /* both subdomains on rank 0 */
82*c4762a1bSJed Brown   }
83*c4762a1bSJed Brown   if (rank) {
84*c4762a1bSJed Brown     jlow = Jlow+1; jhigh = Jhigh+1;
85*c4762a1bSJed Brown   }
86*c4762a1bSJed Brown   else {
87*c4762a1bSJed Brown     jlow = Jlow; jhigh = Jhigh;
88*c4762a1bSJed Brown   }
89*c4762a1bSJed Brown   sort_rows = PETSC_FALSE;
90*c4762a1bSJed Brown   ierr      = PetscOptionsGetBool(NULL,NULL, "-sort_rows", &sort_rows, NULL);CHKERRQ(ierr);
91*c4762a1bSJed Brown   sort_cols = PETSC_FALSE;
92*c4762a1bSJed Brown   ierr      = PetscOptionsGetBool(NULL,NULL, "-sort_cols", &sort_cols, NULL);CHKERRQ(ierr);
93*c4762a1bSJed Brown   for (l = 0; l < nsub; ++l) {
94*c4762a1bSJed Brown     ierr = PetscMalloc1(12, &subindices);CHKERRQ(ierr);
95*c4762a1bSJed Brown     k    = 0;
96*c4762a1bSJed Brown     for (i = 0; i < 4; ++i) {
97*c4762a1bSJed Brown       for (j = jlow[l]; j < jhigh[l]; ++j) {
98*c4762a1bSJed Brown         subindices[k] = j*4+i;
99*c4762a1bSJed Brown         k++;
100*c4762a1bSJed Brown       }
101*c4762a1bSJed Brown     }
102*c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF, 12, subindices, PETSC_OWN_POINTER, rowis+l);CHKERRQ(ierr);
103*c4762a1bSJed Brown     if ((sort_rows && !sort_cols) || (!sort_rows && sort_cols)) {
104*c4762a1bSJed Brown       ierr = ISDuplicate(rowis[l],colis+l);CHKERRQ(ierr);
105*c4762a1bSJed Brown     } else {
106*c4762a1bSJed Brown       ierr = PetscObjectReference((PetscObject)rowis[l]);CHKERRQ(ierr);
107*c4762a1bSJed Brown       colis[l] = rowis[l];
108*c4762a1bSJed Brown     }
109*c4762a1bSJed Brown     if (sort_rows) {
110*c4762a1bSJed Brown       ierr = ISSort(rowis[l]);CHKERRQ(ierr);
111*c4762a1bSJed Brown     }
112*c4762a1bSJed Brown     if (sort_cols) {
113*c4762a1bSJed Brown       ierr = ISSort(colis[l]);CHKERRQ(ierr);
114*c4762a1bSJed Brown     }
115*c4762a1bSJed Brown   }
116*c4762a1bSJed Brown 
117*c4762a1bSJed Brown   ierr = MatCreateSubMatrices(A,nsub,rowis,colis,MAT_INITIAL_MATRIX, &S);CHKERRQ(ierr);
118*c4762a1bSJed Brown 
119*c4762a1bSJed Brown   show_inversions = PETSC_FALSE;
120*c4762a1bSJed Brown 
121*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL, "-show_inversions", &show_inversions, NULL);CHKERRQ(ierr);
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown   inversions = 0;
124*c4762a1bSJed Brown   for (p = 0; p < size; ++p) {
125*c4762a1bSJed Brown     if (p == rank) {
126*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, "[%D:%D]: Number of subdomains: %D:\n", rank, size, nsub);CHKERRQ(ierr);
127*c4762a1bSJed Brown       for (l = 0; l < nsub; ++l) {
128*c4762a1bSJed Brown         PetscInt i0, i1;
129*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "[%D:%D]: Subdomain row IS %D:\n", rank, size, l);CHKERRQ(ierr);
130*c4762a1bSJed Brown         ierr = ISView(rowis[l],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
131*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "[%D:%D]: Subdomain col IS %D:\n", rank, size, l);CHKERRQ(ierr);
132*c4762a1bSJed Brown         ierr = ISView(colis[l],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
133*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "[%D:%D]: Submatrix %D:\n", rank, size, l);CHKERRQ(ierr);
134*c4762a1bSJed Brown         ierr = MatView(S[l],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
135*c4762a1bSJed Brown         if (show_inversions) {
136*c4762a1bSJed Brown           ierr = MatGetOwnershipRange(S[l], &i0,&i1);CHKERRQ(ierr);
137*c4762a1bSJed Brown           for (i = i0; i < i1; ++i) {
138*c4762a1bSJed Brown             ierr = MatGetRow(S[l], i, &ncols, &cols, NULL);CHKERRQ(ierr);
139*c4762a1bSJed Brown             for (j = 1; j < ncols; ++j) {
140*c4762a1bSJed Brown               if (cols[j] < cols[j-1]) {
141*c4762a1bSJed Brown                 ierr = PetscPrintf(PETSC_COMM_SELF, "***Inversion in row %D: col[%D] = %D < %D = col[%D]\n", i, j, cols[j], cols[j-1], j-1);CHKERRQ(ierr);
142*c4762a1bSJed Brown                 inversions++;
143*c4762a1bSJed Brown               }
144*c4762a1bSJed Brown             }
145*c4762a1bSJed Brown             ierr = MatRestoreRow(S[l], i, &ncols, &cols, NULL);CHKERRQ(ierr);
146*c4762a1bSJed Brown           }
147*c4762a1bSJed Brown         }
148*c4762a1bSJed Brown       }
149*c4762a1bSJed Brown     }
150*c4762a1bSJed Brown     ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr);
151*c4762a1bSJed Brown   }
152*c4762a1bSJed Brown   if (show_inversions) {
153*c4762a1bSJed Brown     ierr = MPI_Reduce(&inversions,&total_inversions,1,MPIU_INT, MPI_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
154*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "*Total inversions: %D\n", total_inversions);CHKERRQ(ierr);
155*c4762a1bSJed Brown   }
156*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown   for (l = 0; l < nsub; ++l) {
159*c4762a1bSJed Brown     ierr = ISDestroy(&(rowis[l]));CHKERRQ(ierr);
160*c4762a1bSJed Brown     ierr = ISDestroy(&(colis[l]));CHKERRQ(ierr);
161*c4762a1bSJed Brown   }
162*c4762a1bSJed Brown   ierr = MatDestroySubMatrices(nsub,&S);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = PetscFinalize();
164*c4762a1bSJed Brown   return ierr;
165*c4762a1bSJed Brown }
166