1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 22c0dbf93SJed Brown 32c0dbf93SJed Brown typedef struct { 42c0dbf93SJed Brown Mat Top; 589f730bfSLawrence Mitchell PetscBool rowisblock; 689f730bfSLawrence Mitchell PetscBool colisblock; 72c0dbf93SJed Brown PetscErrorCode (*SetValues)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 82c0dbf93SJed Brown PetscErrorCode (*SetValuesBlocked)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 92c0dbf93SJed Brown } Mat_LocalRef; 102c0dbf93SJed Brown 112c0dbf93SJed Brown /* These need to be macros because they use sizeof */ 129371c9d4SSatish Balay #define IndexSpaceGet(buf, nrow, ncol, irowm, icolm) \ 139371c9d4SSatish Balay do { \ 14dd39110bSPierre Jolivet if (nrow + ncol > (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) { \ 159566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrow, &irowm, ncol, &icolm)); \ 162c0dbf93SJed Brown } else { \ 172c0dbf93SJed Brown irowm = &buf[0]; \ 182c0dbf93SJed Brown icolm = &buf[nrow]; \ 192c0dbf93SJed Brown } \ 202c0dbf93SJed Brown } while (0) 212c0dbf93SJed Brown 229371c9d4SSatish Balay #define IndexSpaceRestore(buf, nrow, ncol, irowm, icolm) \ 239371c9d4SSatish Balay do { \ 2448a46eb9SPierre Jolivet if (nrow + ncol > (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) PetscCall(PetscFree2(irowm, icolm)); \ 252c0dbf93SJed Brown } while (0) 262c0dbf93SJed Brown 27d71ae5a4SJacob Faibussowitsch static void BlockIndicesExpand(PetscInt n, const PetscInt idx[], PetscInt bs, PetscInt idxm[]) 28d71ae5a4SJacob Faibussowitsch { 292c0dbf93SJed Brown PetscInt i, j; 302c0dbf93SJed Brown for (i = 0; i < n; i++) { 31ad540459SPierre Jolivet for (j = 0; j < bs; j++) idxm[i * bs + j] = idx[i] * bs + j; 322c0dbf93SJed Brown } 332c0dbf93SJed Brown } 342c0dbf93SJed Brown 35d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Block(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) 36d71ae5a4SJacob Faibussowitsch { 372c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef *)A->data; 38a497884dSSatish Balay PetscInt buf[4096], *irowm = NULL, *icolm; /* suppress maybe-uninitialized warning */ 392c0dbf93SJed Brown 402c0dbf93SJed Brown PetscFunctionBegin; 413ba16761SJacob Faibussowitsch if (!nrow || !ncol) PetscFunctionReturn(PETSC_SUCCESS); 422c0dbf93SJed Brown IndexSpaceGet(buf, nrow, ncol, irowm, icolm); 439566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow, irow, irowm)); 449566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol, icol, icolm)); 459566063dSJacob Faibussowitsch PetscCall((*lr->SetValuesBlocked)(lr->Top, nrow, irowm, ncol, icolm, y, addv)); 462c0dbf93SJed Brown IndexSpaceRestore(buf, nrow, ncol, irowm, icolm); 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 482c0dbf93SJed Brown } 492c0dbf93SJed Brown 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Scalar(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) 51d71ae5a4SJacob Faibussowitsch { 522c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef *)A->data; 536e191704SLawrence Mitchell PetscInt rbs, cbs, buf[4096], *irowm, *icolm; 542c0dbf93SJed Brown 552c0dbf93SJed Brown PetscFunctionBegin; 569566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, &rbs, &cbs)); 576e191704SLawrence Mitchell IndexSpaceGet(buf, nrow * rbs, ncol * cbs, irowm, icolm); 586e191704SLawrence Mitchell BlockIndicesExpand(nrow, irow, rbs, irowm); 596e191704SLawrence Mitchell BlockIndicesExpand(ncol, icol, cbs, icolm); 609566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow * rbs, irowm, irowm)); 619566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol * cbs, icolm, icolm)); 629566063dSJacob Faibussowitsch PetscCall((*lr->SetValues)(lr->Top, nrow * rbs, irowm, ncol * cbs, icolm, y, addv)); 636e191704SLawrence Mitchell IndexSpaceRestore(buf, nrow * rbs, ncol * cbs, irowm, icolm); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 652c0dbf93SJed Brown } 662c0dbf93SJed Brown 67d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) 68d71ae5a4SJacob Faibussowitsch { 692c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef *)A->data; 702c0dbf93SJed Brown PetscInt buf[4096], *irowm, *icolm; 712c0dbf93SJed Brown 722c0dbf93SJed Brown PetscFunctionBegin; 732c0dbf93SJed Brown IndexSpaceGet(buf, nrow, ncol, irowm, icolm); 7489f730bfSLawrence Mitchell /* If the row IS defining this submatrix was an ISBLOCK, then the unblocked LGMapApply is the right one to use. If 7589f730bfSLawrence Mitchell * instead it was (say) an ISSTRIDE with a block size > 1, then we need to use LGMapApplyBlock */ 7689f730bfSLawrence Mitchell if (lr->rowisblock) { 779566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(A->rmap->mapping, nrow, irow, irowm)); 7889f730bfSLawrence Mitchell } else { 799566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow, irow, irowm)); 8089f730bfSLawrence Mitchell } 8189f730bfSLawrence Mitchell /* As above, but for the column IS. */ 8289f730bfSLawrence Mitchell if (lr->colisblock) { 839566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(A->cmap->mapping, ncol, icol, icolm)); 8489f730bfSLawrence Mitchell } else { 859566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol, icol, icolm)); 8689f730bfSLawrence Mitchell } 879566063dSJacob Faibussowitsch PetscCall((*lr->SetValues)(lr->Top, nrow, irowm, ncol, icolm, y, addv)); 882c0dbf93SJed Brown IndexSpaceRestore(buf, nrow, ncol, irowm, icolm); 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 902c0dbf93SJed Brown } 912c0dbf93SJed Brown 922265a299SJed Brown /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */ 93d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISL2GCompose(IS is, ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *cltog) 94d71ae5a4SJacob Faibussowitsch { 952c0dbf93SJed Brown const PetscInt *idx; 962c0dbf93SJed Brown PetscInt m, *idxm; 9712e05225SLawrence Mitchell PetscInt bs; 98565245c5SBarry Smith PetscBool isblock; 992c0dbf93SJed Brown 1002c0dbf93SJed Brown PetscFunctionBegin; 101b3e8af9fSJed Brown PetscValidHeaderSpecific(is, IS_CLASSID, 1); 102b3e8af9fSJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 2); 1034f572ea9SToby Isaac PetscAssertPointer(cltog, 3); 1049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock)); 105565245c5SBarry Smith if (isblock) { 10612e05225SLawrence Mitchell PetscInt lbs; 107565245c5SBarry Smith 1089566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is, &bs)); 1099566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog, &lbs)); 110565245c5SBarry Smith if (bs == lbs) { 1119566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &m)); 112565245c5SBarry Smith m = m / bs; 1139566063dSJacob Faibussowitsch PetscCall(ISBlockGetIndices(is, &idx)); 1149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idxm)); 1159566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(ltog, m, idx, idxm)); 1169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog)); 1179566063dSJacob Faibussowitsch PetscCall(ISBlockRestoreIndices(is, &idx)); 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119565245c5SBarry Smith } 120565245c5SBarry Smith } 1219566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &m)); 1229566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idx)); 1239566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is, &bs)); 1249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idxm)); 1252c0dbf93SJed Brown if (ltog) { 1269566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltog, m, idx, idxm)); 1272c0dbf93SJed Brown } else { 1289566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(idxm, idx, m)); 1292c0dbf93SJed Brown } 1309566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog)); 1319566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idx)); 1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1332c0dbf93SJed Brown } 1342c0dbf93SJed Brown 135d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISL2GComposeBlock(IS is, ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *cltog) 136d71ae5a4SJacob Faibussowitsch { 1372c0dbf93SJed Brown const PetscInt *idx; 13845b6f7e9SBarry Smith PetscInt m, *idxm, bs; 1392c0dbf93SJed Brown 1402c0dbf93SJed Brown PetscFunctionBegin; 141b3e8af9fSJed Brown PetscValidHeaderSpecific(is, IS_CLASSID, 1); 142b3e8af9fSJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 2); 1434f572ea9SToby Isaac PetscAssertPointer(cltog, 3); 1449566063dSJacob Faibussowitsch PetscCall(ISBlockGetLocalSize(is, &m)); 1459566063dSJacob Faibussowitsch PetscCall(ISBlockGetIndices(is, &idx)); 1469566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog, &bs)); 1479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idxm)); 1482c0dbf93SJed Brown if (ltog) { 1499566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(ltog, m, idx, idxm)); 1502c0dbf93SJed Brown } else { 1519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(idxm, idx, m)); 1522c0dbf93SJed Brown } 1539566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog)); 1549566063dSJacob Faibussowitsch PetscCall(ISBlockRestoreIndices(is, &idx)); 1553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1562c0dbf93SJed Brown } 1572c0dbf93SJed Brown 158*06024a9cSStefano Zampini static PetscErrorCode MatZeroRowsLocal_LocalRef(Mat A, PetscInt n, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 159*06024a9cSStefano Zampini { 160*06024a9cSStefano Zampini PetscInt *rows_l; 161*06024a9cSStefano Zampini Mat_LocalRef *lr = (Mat_LocalRef *)A->data; 162*06024a9cSStefano Zampini 163*06024a9cSStefano Zampini PetscFunctionBegin; 164*06024a9cSStefano Zampini PetscCall(PetscMalloc1(n, &rows_l)); 165*06024a9cSStefano Zampini PetscCall(ISLocalToGlobalMappingApply(A->rmap->mapping, n, rows, rows_l)); 166*06024a9cSStefano Zampini PetscCall(MatZeroRows(lr->Top, n, rows_l, diag, x, b)); 167*06024a9cSStefano Zampini PetscCall(PetscFree(rows_l)); 168*06024a9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 169*06024a9cSStefano Zampini } 170*06024a9cSStefano Zampini 171*06024a9cSStefano Zampini static PetscErrorCode MatZeroRowsColumnsLocal_LocalRef(Mat A, PetscInt n, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 172*06024a9cSStefano Zampini { 173*06024a9cSStefano Zampini PetscInt *rows_l; 174*06024a9cSStefano Zampini Mat_LocalRef *lr = (Mat_LocalRef *)A->data; 175*06024a9cSStefano Zampini 176*06024a9cSStefano Zampini PetscFunctionBegin; 177*06024a9cSStefano Zampini PetscCall(PetscMalloc1(n, &rows_l)); 178*06024a9cSStefano Zampini PetscCall(ISLocalToGlobalMappingApply(A->rmap->mapping, n, rows, rows_l)); 179*06024a9cSStefano Zampini PetscCall(MatZeroRowsColumns(lr->Top, n, rows_l, diag, x, b)); 180*06024a9cSStefano Zampini PetscCall(PetscFree(rows_l)); 181*06024a9cSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 182*06024a9cSStefano Zampini } 183*06024a9cSStefano Zampini 184d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_LocalRef(Mat B) 185d71ae5a4SJacob Faibussowitsch { 1862c0dbf93SJed Brown PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(PetscFree(B->data)); 1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1892c0dbf93SJed Brown } 1902c0dbf93SJed Brown 1912c0dbf93SJed Brown /*@ 19211a5261eSBarry Smith MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly, that is to set values into the matrix 1932c0dbf93SJed Brown 1942c0dbf93SJed Brown Not Collective 1952c0dbf93SJed Brown 1964165533cSJose E. Roman Input Parameters: 1972ef1f0ffSBarry Smith + A - full matrix, generally parallel 1982c0dbf93SJed Brown . isrow - Local index set for the rows 1992c0dbf93SJed Brown - iscol - Local index set for the columns 2002c0dbf93SJed Brown 2014165533cSJose E. Roman Output Parameter: 2022ef1f0ffSBarry Smith . newmat - new serial `Mat` 2032c0dbf93SJed Brown 2042c0dbf93SJed Brown Level: developer 2052c0dbf93SJed Brown 2062c0dbf93SJed Brown Notes: 20711a5261eSBarry Smith Most will use `MatGetLocalSubMatrix()` which returns a real matrix corresponding to the local 2082c0dbf93SJed Brown block if it available, such as with matrix formats that store these blocks separately. 2092c0dbf93SJed Brown 21011a5261eSBarry Smith The new matrix forwards `MatSetValuesLocal()` and `MatSetValuesBlockedLocal()` to the global system. 21111a5261eSBarry Smith In general, it does not define `MatMult()` or any other functions. Local submatrices can be nested. 2122c0dbf93SJed Brown 213fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATSUBMATRIX`, `MatCreateSubMatrixVirtual()`, `MatSetValuesLocal()`, `MatSetValuesBlockedLocal()`, `MatGetLocalSubMatrix()`, `MatCreateSubMatrix()` 2142c0dbf93SJed Brown @*/ 215d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateLocalRef(Mat A, IS isrow, IS iscol, Mat *newmat) 216d71ae5a4SJacob Faibussowitsch { 2172c0dbf93SJed Brown Mat_LocalRef *lr; 2182c0dbf93SJed Brown Mat B; 2192c0dbf93SJed Brown PetscInt m, n; 2202c0dbf93SJed Brown PetscBool islr; 2212c0dbf93SJed Brown 2222c0dbf93SJed Brown PetscFunctionBegin; 2232c0dbf93SJed Brown PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2242c0dbf93SJed Brown PetscValidHeaderSpecific(isrow, IS_CLASSID, 2); 2252c0dbf93SJed Brown PetscValidHeaderSpecific(iscol, IS_CLASSID, 3); 2264f572ea9SToby Isaac PetscAssertPointer(newmat, 4); 22728b400f6SJacob Faibussowitsch PetscCheck(A->rmap->mapping, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Matrix must have local to global mapping provided before this call"); 228f4259b30SLisandro Dalcin *newmat = NULL; 2292c0dbf93SJed Brown 2309566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &B)); 2319566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow, &m)); 2329566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscol, &n)); 2339566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, n, m, n)); 2349566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATLOCALREF)); 2359566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 2362c0dbf93SJed Brown 2372c0dbf93SJed Brown B->ops->destroy = MatDestroy_LocalRef; 2382c0dbf93SJed Brown 2394dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&lr)); 2402c0dbf93SJed Brown B->data = (void *)lr; 2412c0dbf93SJed Brown 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATLOCALREF, &islr)); 2432265a299SJed Brown if (islr) { 2442265a299SJed Brown Mat_LocalRef *alr = (Mat_LocalRef *)A->data; 2452265a299SJed Brown lr->Top = alr->Top; 2462265a299SJed Brown } else { 2472265a299SJed Brown /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */ 2482265a299SJed Brown lr->Top = A; 2492265a299SJed Brown } 2502265a299SJed Brown { 2512c0dbf93SJed Brown ISLocalToGlobalMapping rltog, cltog; 2526e191704SLawrence Mitchell PetscInt arbs, acbs, rbs, cbs; 2532c0dbf93SJed Brown 2542265a299SJed Brown /* We will translate directly to global indices for the top level */ 2552c0dbf93SJed Brown lr->SetValues = MatSetValues; 2562c0dbf93SJed Brown lr->SetValuesBlocked = MatSetValuesBlocked; 2572c0dbf93SJed Brown 2582c0dbf93SJed Brown B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar; 259*06024a9cSStefano Zampini B->ops->zerorowslocal = MatZeroRowsLocal_LocalRef; 260*06024a9cSStefano Zampini B->ops->zerorowscolumnslocal = MatZeroRowsColumnsLocal_LocalRef; 2612205254eSKarl Rupp 2629566063dSJacob Faibussowitsch PetscCall(ISL2GCompose(isrow, A->rmap->mapping, &rltog)); 263992144d0SBarry Smith if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) { 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rltog)); 2652c0dbf93SJed Brown cltog = rltog; 2662c0dbf93SJed Brown } else { 2679566063dSJacob Faibussowitsch PetscCall(ISL2GCompose(iscol, A->cmap->mapping, &cltog)); 2682c0dbf93SJed Brown } 26989f730bfSLawrence Mitchell /* Remember if the ISes we used to pull out the submatrix are of type ISBLOCK. Will be used later in 27089f730bfSLawrence Mitchell * MatSetValuesLocal_LocalRef_Scalar. */ 2719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)isrow, ISBLOCK, &lr->rowisblock)); 2729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISBLOCK, &lr->colisblock)); 2739566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B, rltog, cltog)); 2749566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rltog)); 2759566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cltog)); 2762c0dbf93SJed Brown 2779566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, &arbs, &acbs)); 2789566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(isrow, &rbs)); 2799566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(iscol, &cbs)); 2806e191704SLawrence Mitchell /* Always support block interface insertion on submatrix */ 2819566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(B->rmap, rbs)); 2829566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(B->cmap, cbs)); 2836e191704SLawrence Mitchell if (arbs != rbs || acbs != cbs || (arbs == 1 && acbs == 1)) { 2842c0dbf93SJed Brown /* Top-level matrix has different block size, so we have to call its scalar insertion interface */ 2852c0dbf93SJed Brown B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar; 2862c0dbf93SJed Brown } else { 2872c0dbf93SJed Brown /* Block sizes match so we can forward values to the top level using the block interface */ 2882c0dbf93SJed Brown B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block; 28926fbe8dcSKarl Rupp 2909566063dSJacob Faibussowitsch PetscCall(ISL2GComposeBlock(isrow, A->rmap->mapping, &rltog)); 29145b6f7e9SBarry Smith if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) { 2929566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rltog)); 2932c0dbf93SJed Brown cltog = rltog; 2942c0dbf93SJed Brown } else { 2959566063dSJacob Faibussowitsch PetscCall(ISL2GComposeBlock(iscol, A->cmap->mapping, &cltog)); 2962c0dbf93SJed Brown } 2979566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B, rltog, cltog)); 2989566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rltog)); 2999566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cltog)); 3002c0dbf93SJed Brown } 3012c0dbf93SJed Brown } 3022c0dbf93SJed Brown *newmat = B; 3033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3042c0dbf93SJed Brown } 305