1 2 const char help[] = "Test MatCreateLocalRef()\n\n"; 3 4 #include <petscmat.h> 5 6 static PetscErrorCode GetLocalRef(Mat A, IS isrow, IS iscol, Mat *B) { 7 IS istmp; 8 9 PetscFunctionBegin; 10 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Extracting LocalRef with isrow:\n")); 11 PetscCall(ISOnComm(isrow, PETSC_COMM_WORLD, PETSC_COPY_VALUES, &istmp)); 12 PetscCall(ISView(istmp, PETSC_VIEWER_STDOUT_WORLD)); 13 PetscCall(ISDestroy(&istmp)); 14 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Extracting LocalRef with iscol (only rank=0 shown):\n")); 15 PetscCall(ISOnComm(iscol, PETSC_COMM_WORLD, PETSC_COPY_VALUES, &istmp)); 16 PetscCall(ISView(istmp, PETSC_VIEWER_STDOUT_WORLD)); 17 PetscCall(ISDestroy(&istmp)); 18 19 PetscCall(MatCreateLocalRef(A, isrow, iscol, B)); 20 PetscFunctionReturn(0); 21 } 22 23 int main(int argc, char *argv[]) { 24 MPI_Comm comm; 25 Mat J, B; 26 PetscInt i, j, k, l, rstart, rend, m, n, top_bs, row_bs, col_bs, nlocblocks, *idx, nrowblocks, ncolblocks, *ridx, *cidx, *icol, *irow; 27 PetscScalar *vals; 28 ISLocalToGlobalMapping brmap; 29 IS is0, is1; 30 PetscBool diag, blocked; 31 32 PetscFunctionBeginUser; 33 PetscCall(PetscInitialize(&argc, &argv, 0, help)); 34 comm = PETSC_COMM_WORLD; 35 36 PetscOptionsBegin(comm, NULL, "LocalRef Test Options", NULL); 37 { 38 top_bs = 2; 39 row_bs = 2; 40 col_bs = 2; 41 diag = PETSC_FALSE; 42 blocked = PETSC_FALSE; 43 PetscCall(PetscOptionsInt("-top_bs", "Block size of top-level matrix", 0, top_bs, &top_bs, NULL)); 44 PetscCall(PetscOptionsInt("-row_bs", "Block size of row map", 0, row_bs, &row_bs, NULL)); 45 PetscCall(PetscOptionsInt("-col_bs", "Block size of col map", 0, col_bs, &col_bs, NULL)); 46 PetscCall(PetscOptionsBool("-diag", "Extract a diagonal black", 0, diag, &diag, NULL)); 47 PetscCall(PetscOptionsBool("-blocked", "Use block insertion", 0, blocked, &blocked, NULL)); 48 } 49 PetscOptionsEnd(); 50 51 PetscCall(MatCreate(comm, &J)); 52 PetscCall(MatSetSizes(J, 6, 6, PETSC_DETERMINE, PETSC_DETERMINE)); 53 PetscCall(MatSetBlockSize(J, top_bs)); 54 PetscCall(MatSetFromOptions(J)); 55 PetscCall(MatSeqBAIJSetPreallocation(J, top_bs, PETSC_DECIDE, 0)); 56 PetscCall(MatMPIBAIJSetPreallocation(J, top_bs, PETSC_DECIDE, 0, PETSC_DECIDE, 0)); 57 PetscCall(MatSetUp(J)); 58 PetscCall(MatGetSize(J, &m, &n)); 59 PetscCall(MatGetOwnershipRange(J, &rstart, &rend)); 60 61 nlocblocks = (rend - rstart) / top_bs + 2; 62 PetscCall(PetscMalloc1(nlocblocks, &idx)); 63 for (i = 0; i < nlocblocks; i++) { idx[i] = (rstart / top_bs + i - 1 + m / top_bs) % (m / top_bs); } 64 PetscCall(ISLocalToGlobalMappingCreate(comm, top_bs, nlocblocks, idx, PETSC_OWN_POINTER, &brmap)); 65 PetscCall(PetscPrintf(comm, "Block ISLocalToGlobalMapping:\n")); 66 PetscCall(ISLocalToGlobalMappingView(brmap, PETSC_VIEWER_STDOUT_WORLD)); 67 68 PetscCall(MatSetLocalToGlobalMapping(J, brmap, brmap)); 69 PetscCall(ISLocalToGlobalMappingDestroy(&brmap)); 70 71 /* Create index sets for local submatrix */ 72 nrowblocks = (rend - rstart) / row_bs; 73 ncolblocks = (rend - rstart) / col_bs; 74 PetscCall(PetscMalloc2(nrowblocks, &ridx, ncolblocks, &cidx)); 75 for (i = 0; i < nrowblocks; i++) ridx[i] = i + ((i > nrowblocks / 2) ^ !rstart); 76 for (i = 0; i < ncolblocks; i++) cidx[i] = i + ((i < ncolblocks / 2) ^ !!rstart); 77 PetscCall(ISCreateBlock(PETSC_COMM_SELF, row_bs, nrowblocks, ridx, PETSC_COPY_VALUES, &is0)); 78 PetscCall(ISCreateBlock(PETSC_COMM_SELF, col_bs, ncolblocks, cidx, PETSC_COPY_VALUES, &is1)); 79 PetscCall(PetscFree2(ridx, cidx)); 80 81 if (diag) { 82 PetscCall(ISDestroy(&is1)); 83 PetscCall(PetscObjectReference((PetscObject)is0)); 84 is1 = is0; 85 ncolblocks = nrowblocks; 86 } 87 88 PetscCall(GetLocalRef(J, is0, is1, &B)); 89 90 PetscCall(MatZeroEntries(J)); 91 92 PetscCall(PetscMalloc3(row_bs, &irow, col_bs, &icol, row_bs * col_bs, &vals)); 93 for (i = 0; i < nrowblocks; i++) { 94 for (j = 0; j < ncolblocks; j++) { 95 for (k = 0; k < row_bs; k++) { 96 for (l = 0; l < col_bs; l++) { vals[k * col_bs + l] = i * 100000 + j * 1000 + k * 10 + l; } 97 irow[k] = i * row_bs + k; 98 } 99 for (l = 0; l < col_bs; l++) icol[l] = j * col_bs + l; 100 if (blocked) { 101 PetscCall(MatSetValuesBlockedLocal(B, 1, &i, 1, &j, vals, ADD_VALUES)); 102 } else { 103 PetscCall(MatSetValuesLocal(B, row_bs, irow, col_bs, icol, vals, ADD_VALUES)); 104 } 105 } 106 } 107 PetscCall(PetscFree3(irow, icol, vals)); 108 109 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 110 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 111 112 PetscCall(MatView(J, PETSC_VIEWER_STDOUT_WORLD)); 113 114 PetscCall(ISDestroy(&is0)); 115 PetscCall(ISDestroy(&is1)); 116 PetscCall(MatDestroy(&B)); 117 PetscCall(MatDestroy(&J)); 118 PetscCall(PetscFinalize()); 119 return 0; 120 } 121 122 /*TEST 123 124 test: 125 nsize: 2 126 filter: grep -v "type: mpi" 127 args: -blocked {{0 1}} -mat_type {{aij baij}} 128 129 TEST*/ 130