1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Extract submatrices using unsorted indices. For SEQSBAIJ either sort both rows and columns, or sort none.\n\n"; 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Take a 4x4 grid and form a 5-point stencil graph Laplacian over it. 5c4762a1bSJed Brown Partition the grid into two subdomains by splitting into two in the j-direction (slowest varying). 6c4762a1bSJed Brown Impose an overlap of 1 and order the subdomains with the j-direction varying fastest. 7c4762a1bSJed Brown Extract the subdomain submatrices, one per rank. 8c4762a1bSJed Brown */ 9c4762a1bSJed Brown /* Results: 10c4762a1bSJed Brown Sequential: 11c4762a1bSJed Brown - seqaij: will error out, if rows or columns are unsorted 12c4762a1bSJed Brown - seqbaij: will extract submatrices correctly even for unsorted row or column indices 13c4762a1bSJed Brown - seqsbaij: will extract submatrices correctly even for unsorted row and column indices (both must be sorted or not); 14c4762a1bSJed Brown CANNOT automatically report inversions, because MatGetRow is not available. 15c4762a1bSJed Brown MPI: 16c4762a1bSJed Brown - mpiaij: will error out, if columns are unsorted 17c4762a1bSJed Brown - mpibaij: will error out, if columns are unsorted. 18c4762a1bSJed Brown - mpisbaij: will error out, if columns are unsorted; even with unsorted rows will produce correct submatrices; 19c4762a1bSJed Brown CANNOT automatically report inversions, because MatGetRow is not available. 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown 22c4762a1bSJed Brown #include <petscmat.h> 23c4762a1bSJed Brown #include <petscis.h> 24c4762a1bSJed Brown 25c4762a1bSJed Brown int main(int argc,char **args) 26c4762a1bSJed Brown { 27c4762a1bSJed Brown Mat A, *S; 28c4762a1bSJed Brown IS rowis[2], colis[2]; 29c4762a1bSJed Brown PetscInt n,N,i,j,k,l,nsub,Jlow[2] = {0,1}, *jlow, Jhigh[2] = {3,4}, *jhigh, row, col, *subindices, ncols; 30c4762a1bSJed Brown const PetscInt *cols; 31c4762a1bSJed Brown PetscScalar v; 32c4762a1bSJed Brown PetscMPIInt rank, size, p, inversions, total_inversions; 33c4762a1bSJed Brown PetscBool sort_rows, sort_cols, show_inversions; 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 38*be096a46SBarry 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) { 42c4762a1bSJed Brown n = 8; N = 16; 43c4762a1bSJed Brown } else { 44c4762a1bSJed Brown n = 16; N = 16; 45c4762a1bSJed Brown } 469566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,n,n,N,N)); 479566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 489566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Don't care if the entries are set multiple times by different procs. */ 51c4762a1bSJed Brown for (i=0; i<4; ++i) { 52c4762a1bSJed Brown for (j = 0; j<4; ++j) { 53c4762a1bSJed Brown row = j*4+i; 54c4762a1bSJed Brown v = -1.0; 55c4762a1bSJed Brown if (i>0) { 569566063dSJacob Faibussowitsch col = row-1; PetscCall(MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES)); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown if (i<3) { 599566063dSJacob Faibussowitsch col = row+1; PetscCall(MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES)); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown if (j>0) { 629566063dSJacob Faibussowitsch col = row-4; PetscCall(MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES)); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown if (j<3) { 659566063dSJacob Faibussowitsch col = row+4; PetscCall(MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES)); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown v = 4.0; 689566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&row,&v,INSERT_VALUES)); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown } 719566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 729566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original matrix\n")); 749566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown if (size > 1) { 77c4762a1bSJed Brown nsub = 1; /* one subdomain per rank */ 78c4762a1bSJed Brown } 79c4762a1bSJed Brown else { 80c4762a1bSJed Brown nsub = 2; /* both subdomains on rank 0 */ 81c4762a1bSJed Brown } 82c4762a1bSJed Brown if (rank) { 83c4762a1bSJed Brown jlow = Jlow+1; jhigh = Jhigh+1; 84c4762a1bSJed Brown } 85c4762a1bSJed Brown else { 86c4762a1bSJed Brown jlow = Jlow; jhigh = Jhigh; 87c4762a1bSJed Brown } 88c4762a1bSJed Brown sort_rows = PETSC_FALSE; 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL, "-sort_rows", &sort_rows, NULL)); 90c4762a1bSJed Brown sort_cols = PETSC_FALSE; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL, "-sort_cols", &sort_cols, NULL)); 92c4762a1bSJed Brown for (l = 0; l < nsub; ++l) { 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(12, &subindices)); 94c4762a1bSJed Brown k = 0; 95c4762a1bSJed Brown for (i = 0; i < 4; ++i) { 96c4762a1bSJed Brown for (j = jlow[l]; j < jhigh[l]; ++j) { 97c4762a1bSJed Brown subindices[k] = j*4+i; 98c4762a1bSJed Brown k++; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown } 1019566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 12, subindices, PETSC_OWN_POINTER, rowis+l)); 102c4762a1bSJed Brown if ((sort_rows && !sort_cols) || (!sort_rows && sort_cols)) { 1039566063dSJacob Faibussowitsch PetscCall(ISDuplicate(rowis[l],colis+l)); 104c4762a1bSJed Brown } else { 1059566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rowis[l])); 106c4762a1bSJed Brown colis[l] = rowis[l]; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown if (sort_rows) { 1099566063dSJacob Faibussowitsch PetscCall(ISSort(rowis[l])); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown if (sort_cols) { 1129566063dSJacob Faibussowitsch PetscCall(ISSort(colis[l])); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,nsub,rowis,colis,MAT_INITIAL_MATRIX, &S)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown show_inversions = PETSC_FALSE; 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL, "-show_inversions", &show_inversions, NULL)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown inversions = 0; 123c4762a1bSJed Brown for (p = 0; p < size; ++p) { 124c4762a1bSJed Brown if (p == rank) { 1259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Number of subdomains: %" PetscInt_FMT ":\n", rank, size, nsub)); 126c4762a1bSJed Brown for (l = 0; l < nsub; ++l) { 127c4762a1bSJed Brown PetscInt i0, i1; 1289566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Subdomain row IS %" PetscInt_FMT ":\n", rank, size, l)); 1299566063dSJacob Faibussowitsch PetscCall(ISView(rowis[l],PETSC_VIEWER_STDOUT_SELF)); 1309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Subdomain col IS %" PetscInt_FMT ":\n", rank, size, l)); 1319566063dSJacob Faibussowitsch PetscCall(ISView(colis[l],PETSC_VIEWER_STDOUT_SELF)); 1329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Submatrix %" PetscInt_FMT ":\n", rank, size, l)); 1339566063dSJacob Faibussowitsch PetscCall(MatView(S[l],PETSC_VIEWER_STDOUT_SELF)); 134c4762a1bSJed Brown if (show_inversions) { 1359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(S[l], &i0,&i1)); 136c4762a1bSJed Brown for (i = i0; i < i1; ++i) { 1379566063dSJacob Faibussowitsch PetscCall(MatGetRow(S[l], i, &ncols, &cols, NULL)); 138c4762a1bSJed Brown for (j = 1; j < ncols; ++j) { 139c4762a1bSJed Brown if (cols[j] < cols[j-1]) { 1409566063dSJacob 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)); 141c4762a1bSJed Brown inversions++; 142c4762a1bSJed Brown } 143c4762a1bSJed Brown } 1449566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(S[l], i, &ncols, &cols, NULL)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown } 147c4762a1bSJed Brown } 148c4762a1bSJed Brown } 1499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown if (show_inversions) { 1529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Reduce(&inversions,&total_inversions,1,MPIU_INT, MPI_SUM,0,PETSC_COMM_WORLD)); 1539566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "*Total inversions: %" PetscInt_FMT "\n", total_inversions)); 154c4762a1bSJed Brown } 1559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 156c4762a1bSJed Brown 157c4762a1bSJed Brown for (l = 0; l < nsub; ++l) { 1589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&(rowis[l]))); 1599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&(colis[l]))); 160c4762a1bSJed Brown } 1619566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nsub,&S)); 1629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 163b122ec5aSJacob Faibussowitsch return 0; 164c4762a1bSJed Brown } 165