xref: /petsc/src/mat/tests/ex167.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Extract submatrices using unsorted indices. For SEQSBAIJ either sort both rows and columns, or sort none.\n\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    Take a 4x4 grid and form a 5-point stencil graph Laplacian over it.
4c4762a1bSJed Brown    Partition the grid into two subdomains by splitting into two in the j-direction (slowest varying).
5c4762a1bSJed Brown    Impose an overlap of 1 and order the subdomains with the j-direction varying fastest.
6c4762a1bSJed Brown    Extract the subdomain submatrices, one per rank.
7c4762a1bSJed Brown */
8c4762a1bSJed Brown /* Results:
9c4762a1bSJed Brown     Sequential:
10c4762a1bSJed Brown     - seqaij:   will error out, if rows or columns are unsorted
11c4762a1bSJed Brown     - seqbaij:  will extract submatrices correctly even for unsorted row or column indices
12c4762a1bSJed Brown     - seqsbaij: will extract submatrices correctly even for unsorted row and column indices (both must be sorted or not);
13c4762a1bSJed Brown                 CANNOT automatically report inversions, because MatGetRow is not available.
14c4762a1bSJed Brown     MPI:
15c4762a1bSJed Brown     - mpiaij:   will error out, if columns are unsorted
16c4762a1bSJed Brown     - mpibaij:  will error out, if columns are unsorted.
17c4762a1bSJed Brown     - mpisbaij: will error out, if columns are unsorted; even with unsorted rows will produce correct submatrices;
18c4762a1bSJed Brown                 CANNOT automatically report inversions, because MatGetRow is not available.
19c4762a1bSJed Brown */
20c4762a1bSJed Brown 
21c4762a1bSJed Brown #include <petscmat.h>
22c4762a1bSJed Brown #include <petscis.h>
23c4762a1bSJed Brown 
main(int argc,char ** args)24d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
25d71ae5a4SJacob Faibussowitsch {
26c4762a1bSJed Brown   Mat             A, *S;
27c4762a1bSJed Brown   IS              rowis[2], colis[2];
28c4762a1bSJed Brown   PetscInt        n, N, i, j, k, l, nsub, Jlow[2] = {0, 1}, *jlow, Jhigh[2] = {3, 4}, *jhigh, row, col, *subindices, ncols;
29c4762a1bSJed Brown   const PetscInt *cols;
30c4762a1bSJed Brown   PetscScalar     v;
31c4762a1bSJed Brown   PetscMPIInt     rank, size, p, inversions, total_inversions;
32c4762a1bSJed Brown   PetscBool       sort_rows, sort_cols, show_inversions;
33c4762a1bSJed Brown 
34327415f7SBarry Smith   PetscFunctionBeginUser;
35*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
38be096a46SBarry Smith   PetscCheck(size < 3, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "A uniprocessor or two-processor example only.");
39c4762a1bSJed Brown 
409566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
41c4762a1bSJed Brown   if (size > 1) {
429371c9d4SSatish Balay     n = 8;
439371c9d4SSatish Balay     N = 16;
44c4762a1bSJed Brown   } else {
459371c9d4SSatish Balay     n = 16;
469371c9d4SSatish Balay     N = 16;
47c4762a1bSJed Brown   }
489566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, N, N));
499566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
509566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Don't care if the entries are set multiple times by different procs. */
53c4762a1bSJed Brown   for (i = 0; i < 4; ++i) {
54c4762a1bSJed Brown     for (j = 0; j < 4; ++j) {
55c4762a1bSJed Brown       row = j * 4 + i;
56c4762a1bSJed Brown       v   = -1.0;
57c4762a1bSJed Brown       if (i > 0) {
589371c9d4SSatish Balay         col = row - 1;
599371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &row, 1, &col, &v, INSERT_VALUES));
60c4762a1bSJed Brown       }
61c4762a1bSJed Brown       if (i < 3) {
629371c9d4SSatish Balay         col = row + 1;
639371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &row, 1, &col, &v, INSERT_VALUES));
64c4762a1bSJed Brown       }
65c4762a1bSJed Brown       if (j > 0) {
669371c9d4SSatish Balay         col = row - 4;
679371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &row, 1, &col, &v, INSERT_VALUES));
68c4762a1bSJed Brown       }
69c4762a1bSJed Brown       if (j < 3) {
709371c9d4SSatish Balay         col = row + 4;
719371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &row, 1, &col, &v, INSERT_VALUES));
72c4762a1bSJed Brown       }
73c4762a1bSJed Brown       v = 4.0;
749566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &row, 1, &row, &v, INSERT_VALUES));
75c4762a1bSJed Brown     }
76c4762a1bSJed Brown   }
779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
799566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original matrix\n"));
809566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   if (size > 1) {
83c4762a1bSJed Brown     nsub = 1; /* one subdomain per rank */
849371c9d4SSatish Balay   } else {
85c4762a1bSJed Brown     nsub = 2; /* both subdomains on rank 0 */
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown   if (rank) {
889371c9d4SSatish Balay     jlow  = Jlow + 1;
899371c9d4SSatish Balay     jhigh = Jhigh + 1;
909371c9d4SSatish Balay   } else {
919371c9d4SSatish Balay     jlow  = Jlow;
929371c9d4SSatish Balay     jhigh = Jhigh;
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown   sort_rows = PETSC_FALSE;
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-sort_rows", &sort_rows, NULL));
96c4762a1bSJed Brown   sort_cols = PETSC_FALSE;
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-sort_cols", &sort_cols, NULL));
98c4762a1bSJed Brown   for (l = 0; l < nsub; ++l) {
999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(12, &subindices));
100c4762a1bSJed Brown     k = 0;
101c4762a1bSJed Brown     for (i = 0; i < 4; ++i) {
102c4762a1bSJed Brown       for (j = jlow[l]; j < jhigh[l]; ++j) {
103c4762a1bSJed Brown         subindices[k] = j * 4 + i;
104c4762a1bSJed Brown         k++;
105c4762a1bSJed Brown       }
106c4762a1bSJed Brown     }
1079566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 12, subindices, PETSC_OWN_POINTER, rowis + l));
108c4762a1bSJed Brown     if ((sort_rows && !sort_cols) || (!sort_rows && sort_cols)) {
1099566063dSJacob Faibussowitsch       PetscCall(ISDuplicate(rowis[l], colis + l));
110c4762a1bSJed Brown     } else {
1119566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)rowis[l]));
112c4762a1bSJed Brown       colis[l] = rowis[l];
113c4762a1bSJed Brown     }
11448a46eb9SPierre Jolivet     if (sort_rows) PetscCall(ISSort(rowis[l]));
11548a46eb9SPierre Jolivet     if (sort_cols) PetscCall(ISSort(colis[l]));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, nsub, rowis, colis, MAT_INITIAL_MATRIX, &S));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   show_inversions = PETSC_FALSE;
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-show_inversions", &show_inversions, NULL));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   inversions = 0;
125c4762a1bSJed Brown   for (p = 0; p < size; ++p) {
126c4762a1bSJed Brown     if (p == rank) {
1279566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Number of subdomains: %" PetscInt_FMT ":\n", rank, size, nsub));
128c4762a1bSJed Brown       for (l = 0; l < nsub; ++l) {
129c4762a1bSJed Brown         PetscInt i0, i1;
1309566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Subdomain row IS %" PetscInt_FMT ":\n", rank, size, l));
1319566063dSJacob Faibussowitsch         PetscCall(ISView(rowis[l], PETSC_VIEWER_STDOUT_SELF));
1329566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Subdomain col IS %" PetscInt_FMT ":\n", rank, size, l));
1339566063dSJacob Faibussowitsch         PetscCall(ISView(colis[l], PETSC_VIEWER_STDOUT_SELF));
1349566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Submatrix %" PetscInt_FMT ":\n", rank, size, l));
1359566063dSJacob Faibussowitsch         PetscCall(MatView(S[l], PETSC_VIEWER_STDOUT_SELF));
136c4762a1bSJed Brown         if (show_inversions) {
1379566063dSJacob Faibussowitsch           PetscCall(MatGetOwnershipRange(S[l], &i0, &i1));
138c4762a1bSJed Brown           for (i = i0; i < i1; ++i) {
1399566063dSJacob Faibussowitsch             PetscCall(MatGetRow(S[l], i, &ncols, &cols, NULL));
140c4762a1bSJed Brown             for (j = 1; j < ncols; ++j) {
141c4762a1bSJed Brown               if (cols[j] < cols[j - 1]) {
1429566063dSJacob Faibussowitsch                 PetscCall(PetscPrintf(PETSC_COMM_SELF, "***Inversion in row %" PetscInt_FMT ": col[%" PetscInt_FMT "] = %" PetscInt_FMT " < %" PetscInt_FMT " = col[%" PetscInt_FMT "]\n", i, j, cols[j], cols[j - 1], j - 1));
143c4762a1bSJed Brown                 inversions++;
144c4762a1bSJed Brown               }
145c4762a1bSJed Brown             }
1469566063dSJacob Faibussowitsch             PetscCall(MatRestoreRow(S[l], i, &ncols, &cols, NULL));
147c4762a1bSJed Brown           }
148c4762a1bSJed Brown         }
149c4762a1bSJed Brown       }
150c4762a1bSJed Brown     }
1519566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
152c4762a1bSJed Brown   }
153c4762a1bSJed Brown   if (show_inversions) {
1549566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Reduce(&inversions, &total_inversions, 1, MPIU_INT, MPI_SUM, 0, PETSC_COMM_WORLD));
1559566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "*Total inversions: %" PetscInt_FMT "\n", total_inversions));
156c4762a1bSJed Brown   }
1579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   for (l = 0; l < nsub; ++l) {
160f4f49eeaSPierre Jolivet     PetscCall(ISDestroy(&rowis[l]));
161f4f49eeaSPierre Jolivet     PetscCall(ISDestroy(&colis[l]));
162c4762a1bSJed Brown   }
1639566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(nsub, &S));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
165b122ec5aSJacob Faibussowitsch   return 0;
166c4762a1bSJed Brown }
167