1c4762a1bSJed Brown 2c4762a1bSJed Brown const char help[] = "Test MatCreateLocalRef()\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6*9371c9d4SSatish Balay static PetscErrorCode GetLocalRef(Mat A, IS isrow, IS iscol, Mat *B) { 7c4762a1bSJed Brown IS istmp; 8c4762a1bSJed Brown 9c4762a1bSJed Brown PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Extracting LocalRef with isrow:\n")); 119566063dSJacob Faibussowitsch PetscCall(ISOnComm(isrow, PETSC_COMM_WORLD, PETSC_COPY_VALUES, &istmp)); 129566063dSJacob Faibussowitsch PetscCall(ISView(istmp, PETSC_VIEWER_STDOUT_WORLD)); 139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&istmp)); 149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Extracting LocalRef with iscol (only rank=0 shown):\n")); 159566063dSJacob Faibussowitsch PetscCall(ISOnComm(iscol, PETSC_COMM_WORLD, PETSC_COPY_VALUES, &istmp)); 169566063dSJacob Faibussowitsch PetscCall(ISView(istmp, PETSC_VIEWER_STDOUT_WORLD)); 179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&istmp)); 18c4762a1bSJed Brown 199566063dSJacob Faibussowitsch PetscCall(MatCreateLocalRef(A, isrow, iscol, B)); 20c4762a1bSJed Brown PetscFunctionReturn(0); 21c4762a1bSJed Brown } 22c4762a1bSJed Brown 23*9371c9d4SSatish Balay int main(int argc, char *argv[]) { 24c4762a1bSJed Brown MPI_Comm comm; 25c4762a1bSJed Brown Mat J, B; 26c4762a1bSJed Brown PetscInt i, j, k, l, rstart, rend, m, n, top_bs, row_bs, col_bs, nlocblocks, *idx, nrowblocks, ncolblocks, *ridx, *cidx, *icol, *irow; 27c4762a1bSJed Brown PetscScalar *vals; 28c4762a1bSJed Brown ISLocalToGlobalMapping brmap; 29c4762a1bSJed Brown IS is0, is1; 30c4762a1bSJed Brown PetscBool diag, blocked; 31c4762a1bSJed Brown 32327415f7SBarry Smith PetscFunctionBeginUser; 339566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 34c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 35c4762a1bSJed Brown 36d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "LocalRef Test Options", NULL); 37c4762a1bSJed Brown { 38*9371c9d4SSatish Balay top_bs = 2; 39*9371c9d4SSatish Balay row_bs = 2; 40*9371c9d4SSatish Balay col_bs = 2; 41*9371c9d4SSatish Balay diag = PETSC_FALSE; 42*9371c9d4SSatish Balay blocked = PETSC_FALSE; 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-top_bs", "Block size of top-level matrix", 0, top_bs, &top_bs, NULL)); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-row_bs", "Block size of row map", 0, row_bs, &row_bs, NULL)); 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-col_bs", "Block size of col map", 0, col_bs, &col_bs, NULL)); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-diag", "Extract a diagonal black", 0, diag, &diag, NULL)); 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-blocked", "Use block insertion", 0, blocked, &blocked, NULL)); 48c4762a1bSJed Brown } 49d0609cedSBarry Smith PetscOptionsEnd(); 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &J)); 529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, 6, 6, PETSC_DETERMINE, PETSC_DETERMINE)); 539566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(J, top_bs)); 549566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 559566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(J, top_bs, PETSC_DECIDE, 0)); 569566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(J, top_bs, PETSC_DECIDE, 0, PETSC_DECIDE, 0)); 579566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 589566063dSJacob Faibussowitsch PetscCall(MatGetSize(J, &m, &n)); 599566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(J, &rstart, &rend)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown nlocblocks = (rend - rstart) / top_bs + 2; 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlocblocks, &idx)); 63*9371c9d4SSatish Balay for (i = 0; i < nlocblocks; i++) { idx[i] = (rstart / top_bs + i - 1 + m / top_bs) % (m / top_bs); } 649566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, top_bs, nlocblocks, idx, PETSC_OWN_POINTER, &brmap)); 659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Block ISLocalToGlobalMapping:\n")); 669566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingView(brmap, PETSC_VIEWER_STDOUT_WORLD)); 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(J, brmap, brmap)); 699566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&brmap)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Create index sets for local submatrix */ 72c4762a1bSJed Brown nrowblocks = (rend - rstart) / row_bs; 73c4762a1bSJed Brown ncolblocks = (rend - rstart) / col_bs; 749566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrowblocks, &ridx, ncolblocks, &cidx)); 75c4762a1bSJed Brown for (i = 0; i < nrowblocks; i++) ridx[i] = i + ((i > nrowblocks / 2) ^ !rstart); 76c4762a1bSJed Brown for (i = 0; i < ncolblocks; i++) cidx[i] = i + ((i < ncolblocks / 2) ^ !!rstart); 779566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, row_bs, nrowblocks, ridx, PETSC_COPY_VALUES, &is0)); 789566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF, col_bs, ncolblocks, cidx, PETSC_COPY_VALUES, &is1)); 799566063dSJacob Faibussowitsch PetscCall(PetscFree2(ridx, cidx)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown if (diag) { 829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 839566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)is0)); 84c4762a1bSJed Brown is1 = is0; 85c4762a1bSJed Brown ncolblocks = nrowblocks; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(GetLocalRef(J, is0, is1, &B)); 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J)); 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(row_bs, &irow, col_bs, &icol, row_bs * col_bs, &vals)); 93c4762a1bSJed Brown for (i = 0; i < nrowblocks; i++) { 94c4762a1bSJed Brown for (j = 0; j < ncolblocks; j++) { 95c4762a1bSJed Brown for (k = 0; k < row_bs; k++) { 96*9371c9d4SSatish Balay for (l = 0; l < col_bs; l++) { vals[k * col_bs + l] = i * 100000 + j * 1000 + k * 10 + l; } 97c4762a1bSJed Brown irow[k] = i * row_bs + k; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown for (l = 0; l < col_bs; l++) icol[l] = j * col_bs + l; 100c4762a1bSJed Brown if (blocked) { 1019566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(B, 1, &i, 1, &j, vals, ADD_VALUES)); 102c4762a1bSJed Brown } else { 1039566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(B, row_bs, irow, col_bs, icol, vals, ADD_VALUES)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown } 106c4762a1bSJed Brown } 1079566063dSJacob Faibussowitsch PetscCall(PetscFree3(irow, icol, vals)); 108c4762a1bSJed Brown 1099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 111c4762a1bSJed Brown 1129566063dSJacob Faibussowitsch PetscCall(MatView(J, PETSC_VIEWER_STDOUT_WORLD)); 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is0)); 1159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 1169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1189566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 119b122ec5aSJacob Faibussowitsch return 0; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown /*TEST 123c4762a1bSJed Brown 124c4762a1bSJed Brown test: 125c4762a1bSJed Brown nsize: 2 126c4762a1bSJed Brown filter: grep -v "type: mpi" 127c4762a1bSJed Brown args: -blocked {{0 1}} -mat_type {{aij baij}} 128c4762a1bSJed Brown 129c4762a1bSJed Brown TEST*/ 130