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 PetscErrorCode ierr; 35c4762a1bSJed Brown 36c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 37ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 38ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 39*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size>2,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG, "A uniprocessor or two-processor example only."); 40c4762a1bSJed Brown 41c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 42c4762a1bSJed Brown if (size > 1) { 43c4762a1bSJed Brown n = 8; N = 16; 44c4762a1bSJed Brown } else { 45c4762a1bSJed Brown n = 16; N = 16; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Don't care if the entries are set multiple times by different procs. */ 52c4762a1bSJed Brown for (i=0; i<4; ++i) { 53c4762a1bSJed Brown for (j = 0; j<4; ++j) { 54c4762a1bSJed Brown row = j*4+i; 55c4762a1bSJed Brown v = -1.0; 56c4762a1bSJed Brown if (i>0) { 57c4762a1bSJed Brown col = row-1; ierr = MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown if (i<3) { 60c4762a1bSJed Brown col = row+1; ierr = MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown if (j>0) { 63c4762a1bSJed Brown col = row-4; ierr = MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown if (j<3) { 66c4762a1bSJed Brown col = row+4; ierr = MatSetValues(A,1,&row,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown v = 4.0; 69c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&row,&v,INSERT_VALUES);CHKERRQ(ierr); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown } 72c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Original matrix\n");CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 76c4762a1bSJed Brown 77c4762a1bSJed Brown if (size > 1) { 78c4762a1bSJed Brown nsub = 1; /* one subdomain per rank */ 79c4762a1bSJed Brown } 80c4762a1bSJed Brown else { 81c4762a1bSJed Brown nsub = 2; /* both subdomains on rank 0 */ 82c4762a1bSJed Brown } 83c4762a1bSJed Brown if (rank) { 84c4762a1bSJed Brown jlow = Jlow+1; jhigh = Jhigh+1; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown else { 87c4762a1bSJed Brown jlow = Jlow; jhigh = Jhigh; 88c4762a1bSJed Brown } 89c4762a1bSJed Brown sort_rows = PETSC_FALSE; 90c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL, "-sort_rows", &sort_rows, NULL);CHKERRQ(ierr); 91c4762a1bSJed Brown sort_cols = PETSC_FALSE; 92c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL, "-sort_cols", &sort_cols, NULL);CHKERRQ(ierr); 93c4762a1bSJed Brown for (l = 0; l < nsub; ++l) { 94c4762a1bSJed Brown ierr = PetscMalloc1(12, &subindices);CHKERRQ(ierr); 95c4762a1bSJed Brown k = 0; 96c4762a1bSJed Brown for (i = 0; i < 4; ++i) { 97c4762a1bSJed Brown for (j = jlow[l]; j < jhigh[l]; ++j) { 98c4762a1bSJed Brown subindices[k] = j*4+i; 99c4762a1bSJed Brown k++; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown } 102c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF, 12, subindices, PETSC_OWN_POINTER, rowis+l);CHKERRQ(ierr); 103c4762a1bSJed Brown if ((sort_rows && !sort_cols) || (!sort_rows && sort_cols)) { 104c4762a1bSJed Brown ierr = ISDuplicate(rowis[l],colis+l);CHKERRQ(ierr); 105c4762a1bSJed Brown } else { 106c4762a1bSJed Brown ierr = PetscObjectReference((PetscObject)rowis[l]);CHKERRQ(ierr); 107c4762a1bSJed Brown colis[l] = rowis[l]; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown if (sort_rows) { 110c4762a1bSJed Brown ierr = ISSort(rowis[l]);CHKERRQ(ierr); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown if (sort_cols) { 113c4762a1bSJed Brown ierr = ISSort(colis[l]);CHKERRQ(ierr); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,nsub,rowis,colis,MAT_INITIAL_MATRIX, &S);CHKERRQ(ierr); 118c4762a1bSJed Brown 119c4762a1bSJed Brown show_inversions = PETSC_FALSE; 120c4762a1bSJed Brown 121c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL, "-show_inversions", &show_inversions, NULL);CHKERRQ(ierr); 122c4762a1bSJed Brown 123c4762a1bSJed Brown inversions = 0; 124c4762a1bSJed Brown for (p = 0; p < size; ++p) { 125c4762a1bSJed Brown if (p == rank) { 1268fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Number of subdomains: %" PetscInt_FMT ":\n", rank, size, nsub);CHKERRQ(ierr); 127c4762a1bSJed Brown for (l = 0; l < nsub; ++l) { 128c4762a1bSJed Brown PetscInt i0, i1; 1298fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Subdomain row IS %" PetscInt_FMT ":\n", rank, size, l);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = ISView(rowis[l],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 1318fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Subdomain col IS %" PetscInt_FMT ":\n", rank, size, l);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = ISView(colis[l],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 1338fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF, "[%" PetscInt_FMT ":%" PetscInt_FMT "]: Submatrix %" PetscInt_FMT ":\n", rank, size, l);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = MatView(S[l],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 135c4762a1bSJed Brown if (show_inversions) { 136c4762a1bSJed Brown ierr = MatGetOwnershipRange(S[l], &i0,&i1);CHKERRQ(ierr); 137c4762a1bSJed Brown for (i = i0; i < i1; ++i) { 138c4762a1bSJed Brown ierr = MatGetRow(S[l], i, &ncols, &cols, NULL);CHKERRQ(ierr); 139c4762a1bSJed Brown for (j = 1; j < ncols; ++j) { 140c4762a1bSJed Brown if (cols[j] < cols[j-1]) { 1418fffc762SJacob Faibussowitsch ierr = 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);CHKERRQ(ierr); 142c4762a1bSJed Brown inversions++; 143c4762a1bSJed Brown } 144c4762a1bSJed Brown } 145c4762a1bSJed Brown ierr = MatRestoreRow(S[l], i, &ncols, &cols, NULL);CHKERRQ(ierr); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown } 148c4762a1bSJed Brown } 149c4762a1bSJed Brown } 150ffc4695bSBarry Smith ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRMPI(ierr); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown if (show_inversions) { 153ffc4695bSBarry Smith ierr = MPI_Reduce(&inversions,&total_inversions,1,MPIU_INT, MPI_SUM,0,PETSC_COMM_WORLD);CHKERRMPI(ierr); 1548fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD, "*Total inversions: %" PetscInt_FMT "\n", total_inversions);CHKERRQ(ierr); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 157c4762a1bSJed Brown 158c4762a1bSJed Brown for (l = 0; l < nsub; ++l) { 159c4762a1bSJed Brown ierr = ISDestroy(&(rowis[l]));CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = ISDestroy(&(colis[l]));CHKERRQ(ierr); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown ierr = MatDestroySubMatrices(nsub,&S);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = PetscFinalize(); 164c4762a1bSJed Brown return ierr; 165c4762a1bSJed Brown } 166