1c4762a1bSJed Brown 2c4762a1bSJed Brown const char help[] = "Test MatCreateLocalRef()\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode GetLocalRef(Mat A, IS isrow, IS iscol, Mat *B) 7d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown IS istmp; 9c4762a1bSJed Brown 10c4762a1bSJed Brown PetscFunctionBegin; 119566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Extracting LocalRef with isrow:\n")); 129566063dSJacob Faibussowitsch PetscCall(ISOnComm(isrow, PETSC_COMM_WORLD, PETSC_COPY_VALUES, &istmp)); 139566063dSJacob Faibussowitsch PetscCall(ISView(istmp, PETSC_VIEWER_STDOUT_WORLD)); 149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&istmp)); 159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Extracting LocalRef with iscol (only rank=0 shown):\n")); 169566063dSJacob Faibussowitsch PetscCall(ISOnComm(iscol, PETSC_COMM_WORLD, PETSC_COPY_VALUES, &istmp)); 179566063dSJacob Faibussowitsch PetscCall(ISView(istmp, PETSC_VIEWER_STDOUT_WORLD)); 189566063dSJacob Faibussowitsch PetscCall(ISDestroy(&istmp)); 19c4762a1bSJed Brown 209566063dSJacob Faibussowitsch PetscCall(MatCreateLocalRef(A, isrow, iscol, B)); 21*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22c4762a1bSJed Brown } 23c4762a1bSJed Brown 24d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 25d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown MPI_Comm comm; 27c4762a1bSJed Brown Mat J, B; 28c4762a1bSJed Brown PetscInt i, j, k, l, rstart, rend, m, n, top_bs, row_bs, col_bs, nlocblocks, *idx, nrowblocks, ncolblocks, *ridx, *cidx, *icol, *irow; 29c4762a1bSJed Brown PetscScalar *vals; 30c4762a1bSJed Brown ISLocalToGlobalMapping brmap; 31c4762a1bSJed Brown IS is0, is1; 32c4762a1bSJed Brown PetscBool diag, blocked; 33c4762a1bSJed Brown 34327415f7SBarry Smith PetscFunctionBeginUser; 359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 36c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 37c4762a1bSJed Brown 38d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "LocalRef Test Options", NULL); 39c4762a1bSJed Brown { 409371c9d4SSatish Balay top_bs = 2; 419371c9d4SSatish Balay row_bs = 2; 429371c9d4SSatish Balay col_bs = 2; 439371c9d4SSatish Balay diag = PETSC_FALSE; 449371c9d4SSatish Balay blocked = PETSC_FALSE; 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-top_bs", "Block size of top-level matrix", 0, top_bs, &top_bs, NULL)); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-row_bs", "Block size of row map", 0, row_bs, &row_bs, NULL)); 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-col_bs", "Block size of col map", 0, col_bs, &col_bs, NULL)); 489566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-diag", "Extract a diagonal black", 0, diag, &diag, NULL)); 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-blocked", "Use block insertion", 0, blocked, &blocked, NULL)); 50c4762a1bSJed Brown } 51d0609cedSBarry Smith PetscOptionsEnd(); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &J)); 549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, 6, 6, PETSC_DETERMINE, PETSC_DETERMINE)); 559566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(J, top_bs)); 569566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 579566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(J, top_bs, PETSC_DECIDE, 0)); 589566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(J, top_bs, PETSC_DECIDE, 0, PETSC_DECIDE, 0)); 599566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 609566063dSJacob Faibussowitsch PetscCall(MatGetSize(J, &m, &n)); 619566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(J, &rstart, &rend)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown nlocblocks = (rend - rstart) / top_bs + 2; 649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlocblocks, &idx)); 65ad540459SPierre Jolivet for (i = 0; i < nlocblocks; i++) idx[i] = (rstart / top_bs + i - 1 + m / top_bs) % (m / top_bs); 669566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, top_bs, nlocblocks, idx, PETSC_OWN_POINTER, &brmap)); 679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Block ISLocalToGlobalMapping:\n")); 689566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingView(brmap, PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(J, brmap, brmap)); 719566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&brmap)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* Create index sets for local submatrix */ 74c4762a1bSJed Brown nrowblocks = (rend - rstart) / row_bs; 75c4762a1bSJed Brown ncolblocks = (rend - rstart) / col_bs; 769566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrowblocks, &ridx, ncolblocks, &cidx)); 77c4762a1bSJed Brown for (i = 0; i < nrowblocks; i++) ridx[i] = i + ((i > nrowblocks / 2) ^ !rstart); 78c4762a1bSJed Brown for (i = 0; i < ncolblocks; i++) cidx[i] = i + ((i < ncolblocks / 2) ^ !!rstart); 799566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, row_bs, nrowblocks, ridx, PETSC_COPY_VALUES, &is0)); 809566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, col_bs, ncolblocks, cidx, PETSC_COPY_VALUES, &is1)); 819566063dSJacob Faibussowitsch PetscCall(PetscFree2(ridx, cidx)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown if (diag) { 849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 859566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)is0)); 86c4762a1bSJed Brown is1 = is0; 87c4762a1bSJed Brown ncolblocks = nrowblocks; 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(GetLocalRef(J, is0, is1, &B)); 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J)); 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(row_bs, &irow, col_bs, &icol, row_bs * col_bs, &vals)); 95c4762a1bSJed Brown for (i = 0; i < nrowblocks; i++) { 96c4762a1bSJed Brown for (j = 0; j < ncolblocks; j++) { 97c4762a1bSJed Brown for (k = 0; k < row_bs; k++) { 98ad540459SPierre Jolivet for (l = 0; l < col_bs; l++) vals[k * col_bs + l] = i * 100000 + j * 1000 + k * 10 + l; 99c4762a1bSJed Brown irow[k] = i * row_bs + k; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown for (l = 0; l < col_bs; l++) icol[l] = j * col_bs + l; 102c4762a1bSJed Brown if (blocked) { 1039566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(B, 1, &i, 1, &j, vals, ADD_VALUES)); 104c4762a1bSJed Brown } else { 1059566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(B, row_bs, irow, col_bs, icol, vals, ADD_VALUES)); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 1099566063dSJacob Faibussowitsch PetscCall(PetscFree3(irow, icol, vals)); 110c4762a1bSJed Brown 1119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(MatView(J, PETSC_VIEWER_STDOUT_WORLD)); 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is0)); 1179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 1189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1209566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 121b122ec5aSJacob Faibussowitsch return 0; 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124c4762a1bSJed Brown /*TEST 125c4762a1bSJed Brown 126c4762a1bSJed Brown test: 127c4762a1bSJed Brown nsize: 2 128c4762a1bSJed Brown filter: grep -v "type: mpi" 129c4762a1bSJed Brown args: -blocked {{0 1}} -mat_type {{aij baij}} 130c4762a1bSJed Brown 131c4762a1bSJed Brown TEST*/ 132