12c0dbf93SJed Brown 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 32c0dbf93SJed Brown 42c0dbf93SJed Brown typedef struct { 52c0dbf93SJed Brown Mat Top; 689f730bfSLawrence Mitchell PetscBool rowisblock; 789f730bfSLawrence Mitchell PetscBool colisblock; 82c0dbf93SJed Brown PetscErrorCode (*SetValues)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 92c0dbf93SJed Brown PetscErrorCode (*SetValuesBlocked)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 102c0dbf93SJed Brown } Mat_LocalRef; 112c0dbf93SJed Brown 122c0dbf93SJed Brown /* These need to be macros because they use sizeof */ 139371c9d4SSatish Balay #define IndexSpaceGet(buf, nrow, ncol, irowm, icolm) \ 149371c9d4SSatish Balay do { \ 15dd39110bSPierre Jolivet if (nrow + ncol > (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) { \ 169566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrow, &irowm, ncol, &icolm)); \ 172c0dbf93SJed Brown } else { \ 182c0dbf93SJed Brown irowm = &buf[0]; \ 192c0dbf93SJed Brown icolm = &buf[nrow]; \ 202c0dbf93SJed Brown } \ 212c0dbf93SJed Brown } while (0) 222c0dbf93SJed Brown 239371c9d4SSatish Balay #define IndexSpaceRestore(buf, nrow, ncol, irowm, icolm) \ 249371c9d4SSatish Balay do { \ 2548a46eb9SPierre Jolivet if (nrow + ncol > (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) PetscCall(PetscFree2(irowm, icolm)); \ 262c0dbf93SJed Brown } while (0) 272c0dbf93SJed Brown 28d71ae5a4SJacob Faibussowitsch static void BlockIndicesExpand(PetscInt n, const PetscInt idx[], PetscInt bs, PetscInt idxm[]) 29d71ae5a4SJacob Faibussowitsch { 302c0dbf93SJed Brown PetscInt i, j; 312c0dbf93SJed Brown for (i = 0; i < n; i++) { 32ad540459SPierre Jolivet for (j = 0; j < bs; j++) idxm[i * bs + j] = idx[i] * bs + j; 332c0dbf93SJed Brown } 342c0dbf93SJed Brown } 352c0dbf93SJed Brown 36d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Block(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) 37d71ae5a4SJacob Faibussowitsch { 382c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef *)A->data; 39a497884dSSatish Balay PetscInt buf[4096], *irowm = NULL, *icolm; /* suppress maybe-uninitialized warning */ 402c0dbf93SJed Brown 412c0dbf93SJed Brown PetscFunctionBegin; 42*3ba16761SJacob Faibussowitsch if (!nrow || !ncol) PetscFunctionReturn(PETSC_SUCCESS); 432c0dbf93SJed Brown IndexSpaceGet(buf, nrow, ncol, irowm, icolm); 449566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow, irow, irowm)); 459566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol, icol, icolm)); 469566063dSJacob Faibussowitsch PetscCall((*lr->SetValuesBlocked)(lr->Top, nrow, irowm, ncol, icolm, y, addv)); 472c0dbf93SJed Brown IndexSpaceRestore(buf, nrow, ncol, irowm, icolm); 48*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 492c0dbf93SJed Brown } 502c0dbf93SJed Brown 51d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Scalar(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) 52d71ae5a4SJacob Faibussowitsch { 532c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef *)A->data; 546e191704SLawrence Mitchell PetscInt rbs, cbs, buf[4096], *irowm, *icolm; 552c0dbf93SJed Brown 562c0dbf93SJed Brown PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, &rbs, &cbs)); 586e191704SLawrence Mitchell IndexSpaceGet(buf, nrow * rbs, ncol * cbs, irowm, icolm); 596e191704SLawrence Mitchell BlockIndicesExpand(nrow, irow, rbs, irowm); 606e191704SLawrence Mitchell BlockIndicesExpand(ncol, icol, cbs, icolm); 619566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow * rbs, irowm, irowm)); 629566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol * cbs, icolm, icolm)); 639566063dSJacob Faibussowitsch PetscCall((*lr->SetValues)(lr->Top, nrow * rbs, irowm, ncol * cbs, icolm, y, addv)); 646e191704SLawrence Mitchell IndexSpaceRestore(buf, nrow * rbs, ncol * cbs, irowm, icolm); 65*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 662c0dbf93SJed Brown } 672c0dbf93SJed Brown 68d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) 69d71ae5a4SJacob Faibussowitsch { 702c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef *)A->data; 712c0dbf93SJed Brown PetscInt buf[4096], *irowm, *icolm; 722c0dbf93SJed Brown 732c0dbf93SJed Brown PetscFunctionBegin; 742c0dbf93SJed Brown IndexSpaceGet(buf, nrow, ncol, irowm, icolm); 7589f730bfSLawrence Mitchell /* If the row IS defining this submatrix was an ISBLOCK, then the unblocked LGMapApply is the right one to use. If 7689f730bfSLawrence Mitchell * instead it was (say) an ISSTRIDE with a block size > 1, then we need to use LGMapApplyBlock */ 7789f730bfSLawrence Mitchell if (lr->rowisblock) { 789566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(A->rmap->mapping, nrow, irow, irowm)); 7989f730bfSLawrence Mitchell } else { 809566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow, irow, irowm)); 8189f730bfSLawrence Mitchell } 8289f730bfSLawrence Mitchell /* As above, but for the column IS. */ 8389f730bfSLawrence Mitchell if (lr->colisblock) { 849566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(A->cmap->mapping, ncol, icol, icolm)); 8589f730bfSLawrence Mitchell } else { 869566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol, icol, icolm)); 8789f730bfSLawrence Mitchell } 889566063dSJacob Faibussowitsch PetscCall((*lr->SetValues)(lr->Top, nrow, irowm, ncol, icolm, y, addv)); 892c0dbf93SJed Brown IndexSpaceRestore(buf, nrow, ncol, irowm, icolm); 90*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 912c0dbf93SJed Brown } 922c0dbf93SJed Brown 932265a299SJed Brown /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */ 94d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISL2GCompose(IS is, ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *cltog) 95d71ae5a4SJacob Faibussowitsch { 962c0dbf93SJed Brown const PetscInt *idx; 972c0dbf93SJed Brown PetscInt m, *idxm; 9812e05225SLawrence Mitchell PetscInt bs; 99565245c5SBarry Smith PetscBool isblock; 1002c0dbf93SJed Brown 1012c0dbf93SJed Brown PetscFunctionBegin; 102b3e8af9fSJed Brown PetscValidHeaderSpecific(is, IS_CLASSID, 1); 103b3e8af9fSJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 2); 104b3e8af9fSJed Brown PetscValidPointer(cltog, 3); 1059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock)); 106565245c5SBarry Smith if (isblock) { 10712e05225SLawrence Mitchell PetscInt lbs; 108565245c5SBarry Smith 1099566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is, &bs)); 1109566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog, &lbs)); 111565245c5SBarry Smith if (bs == lbs) { 1129566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &m)); 113565245c5SBarry Smith m = m / bs; 1149566063dSJacob Faibussowitsch PetscCall(ISBlockGetIndices(is, &idx)); 1159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idxm)); 1169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(ltog, m, idx, idxm)); 1179566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog)); 1189566063dSJacob Faibussowitsch PetscCall(ISBlockRestoreIndices(is, &idx)); 119*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120565245c5SBarry Smith } 121565245c5SBarry Smith } 1229566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &m)); 1239566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &idx)); 1249566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(is, &bs)); 1259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idxm)); 1262c0dbf93SJed Brown if (ltog) { 1279566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApply(ltog, m, idx, idxm)); 1282c0dbf93SJed Brown } else { 1299566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(idxm, idx, m)); 1302c0dbf93SJed Brown } 1319566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog)); 1329566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &idx)); 133*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1342c0dbf93SJed Brown } 1352c0dbf93SJed Brown 136d71ae5a4SJacob Faibussowitsch static PetscErrorCode ISL2GComposeBlock(IS is, ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *cltog) 137d71ae5a4SJacob Faibussowitsch { 1382c0dbf93SJed Brown const PetscInt *idx; 13945b6f7e9SBarry Smith PetscInt m, *idxm, bs; 1402c0dbf93SJed Brown 1412c0dbf93SJed Brown PetscFunctionBegin; 142b3e8af9fSJed Brown PetscValidHeaderSpecific(is, IS_CLASSID, 1); 143b3e8af9fSJed Brown PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 2); 144b3e8af9fSJed Brown PetscValidPointer(cltog, 3); 1459566063dSJacob Faibussowitsch PetscCall(ISBlockGetLocalSize(is, &m)); 1469566063dSJacob Faibussowitsch PetscCall(ISBlockGetIndices(is, &idx)); 1479566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog, &bs)); 1489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idxm)); 1492c0dbf93SJed Brown if (ltog) { 1509566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyBlock(ltog, m, idx, idxm)); 1512c0dbf93SJed Brown } else { 1529566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(idxm, idx, m)); 1532c0dbf93SJed Brown } 1549566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog)); 1559566063dSJacob Faibussowitsch PetscCall(ISBlockRestoreIndices(is, &idx)); 156*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1572c0dbf93SJed Brown } 1582c0dbf93SJed Brown 159d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_LocalRef(Mat B) 160d71ae5a4SJacob Faibussowitsch { 1612c0dbf93SJed Brown PetscFunctionBegin; 1629566063dSJacob Faibussowitsch PetscCall(PetscFree(B->data)); 163*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1642c0dbf93SJed Brown } 1652c0dbf93SJed Brown 1662c0dbf93SJed Brown /*@ 16711a5261eSBarry Smith MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly, that is to set values into the matrix 1682c0dbf93SJed Brown 1692c0dbf93SJed Brown Not Collective 1702c0dbf93SJed Brown 1714165533cSJose E. Roman Input Parameters: 1722c0dbf93SJed Brown + A - Full matrix, generally parallel 1732c0dbf93SJed Brown . isrow - Local index set for the rows 1742c0dbf93SJed Brown - iscol - Local index set for the columns 1752c0dbf93SJed Brown 1764165533cSJose E. Roman Output Parameter: 1772c0dbf93SJed Brown . newmat - New serial Mat 1782c0dbf93SJed Brown 1792c0dbf93SJed Brown Level: developer 1802c0dbf93SJed Brown 1812c0dbf93SJed Brown Notes: 18211a5261eSBarry Smith Most will use `MatGetLocalSubMatrix()` which returns a real matrix corresponding to the local 1832c0dbf93SJed Brown block if it available, such as with matrix formats that store these blocks separately. 1842c0dbf93SJed Brown 18511a5261eSBarry Smith The new matrix forwards `MatSetValuesLocal()` and `MatSetValuesBlockedLocal()` to the global system. 18611a5261eSBarry Smith In general, it does not define `MatMult()` or any other functions. Local submatrices can be nested. 1872c0dbf93SJed Brown 18811a5261eSBarry Smith .seealso: MATSUBMATRIX`, `MatCreateSubMatrixVirtual()`, `MatSetValuesLocal()`, `MatSetValuesBlockedLocal()`, `MatGetLocalSubMatrix()`, `MatCreateSubMatrix()` 1892c0dbf93SJed Brown @*/ 190d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateLocalRef(Mat A, IS isrow, IS iscol, Mat *newmat) 191d71ae5a4SJacob Faibussowitsch { 1922c0dbf93SJed Brown Mat_LocalRef *lr; 1932c0dbf93SJed Brown Mat B; 1942c0dbf93SJed Brown PetscInt m, n; 1952c0dbf93SJed Brown PetscBool islr; 1962c0dbf93SJed Brown 1972c0dbf93SJed Brown PetscFunctionBegin; 1982c0dbf93SJed Brown PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 1992c0dbf93SJed Brown PetscValidHeaderSpecific(isrow, IS_CLASSID, 2); 2002c0dbf93SJed Brown PetscValidHeaderSpecific(iscol, IS_CLASSID, 3); 2012c0dbf93SJed Brown PetscValidPointer(newmat, 4); 20228b400f6SJacob Faibussowitsch PetscCheck(A->rmap->mapping, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Matrix must have local to global mapping provided before this call"); 203f4259b30SLisandro Dalcin *newmat = NULL; 2042c0dbf93SJed Brown 2059566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &B)); 2069566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow, &m)); 2079566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscol, &n)); 2089566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, n, m, n)); 2099566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATLOCALREF)); 2109566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 2112c0dbf93SJed Brown 2122c0dbf93SJed Brown B->ops->destroy = MatDestroy_LocalRef; 2132c0dbf93SJed Brown 2144dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&lr)); 2152c0dbf93SJed Brown B->data = (void *)lr; 2162c0dbf93SJed Brown 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATLOCALREF, &islr)); 2182265a299SJed Brown if (islr) { 2192265a299SJed Brown Mat_LocalRef *alr = (Mat_LocalRef *)A->data; 2202265a299SJed Brown lr->Top = alr->Top; 2212265a299SJed Brown } else { 2222265a299SJed Brown /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */ 2232265a299SJed Brown lr->Top = A; 2242265a299SJed Brown } 2252265a299SJed Brown { 2262c0dbf93SJed Brown ISLocalToGlobalMapping rltog, cltog; 2276e191704SLawrence Mitchell PetscInt arbs, acbs, rbs, cbs; 2282c0dbf93SJed Brown 2292265a299SJed Brown /* We will translate directly to global indices for the top level */ 2302c0dbf93SJed Brown lr->SetValues = MatSetValues; 2312c0dbf93SJed Brown lr->SetValuesBlocked = MatSetValuesBlocked; 2322c0dbf93SJed Brown 2332c0dbf93SJed Brown B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar; 2342205254eSKarl Rupp 2359566063dSJacob Faibussowitsch PetscCall(ISL2GCompose(isrow, A->rmap->mapping, &rltog)); 236992144d0SBarry Smith if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) { 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rltog)); 2382c0dbf93SJed Brown cltog = rltog; 2392c0dbf93SJed Brown } else { 2409566063dSJacob Faibussowitsch PetscCall(ISL2GCompose(iscol, A->cmap->mapping, &cltog)); 2412c0dbf93SJed Brown } 24289f730bfSLawrence Mitchell /* Remember if the ISes we used to pull out the submatrix are of type ISBLOCK. Will be used later in 24389f730bfSLawrence Mitchell * MatSetValuesLocal_LocalRef_Scalar. */ 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)isrow, ISBLOCK, &lr->rowisblock)); 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISBLOCK, &lr->colisblock)); 2469566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B, rltog, cltog)); 2479566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rltog)); 2489566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cltog)); 2492c0dbf93SJed Brown 2509566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A, &arbs, &acbs)); 2519566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(isrow, &rbs)); 2529566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(iscol, &cbs)); 2536e191704SLawrence Mitchell /* Always support block interface insertion on submatrix */ 2549566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(B->rmap, rbs)); 2559566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(B->cmap, cbs)); 2566e191704SLawrence Mitchell if (arbs != rbs || acbs != cbs || (arbs == 1 && acbs == 1)) { 2572c0dbf93SJed Brown /* Top-level matrix has different block size, so we have to call its scalar insertion interface */ 2582c0dbf93SJed Brown B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar; 2592c0dbf93SJed Brown } else { 2602c0dbf93SJed Brown /* Block sizes match so we can forward values to the top level using the block interface */ 2612c0dbf93SJed Brown B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block; 26226fbe8dcSKarl Rupp 2639566063dSJacob Faibussowitsch PetscCall(ISL2GComposeBlock(isrow, A->rmap->mapping, &rltog)); 26445b6f7e9SBarry Smith if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) { 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rltog)); 2662c0dbf93SJed Brown cltog = rltog; 2672c0dbf93SJed Brown } else { 2689566063dSJacob Faibussowitsch PetscCall(ISL2GComposeBlock(iscol, A->cmap->mapping, &cltog)); 2692c0dbf93SJed Brown } 2709566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B, rltog, cltog)); 2719566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rltog)); 2729566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cltog)); 2732c0dbf93SJed Brown } 2742c0dbf93SJed Brown } 2752c0dbf93SJed Brown *newmat = B; 276*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2772c0dbf93SJed Brown } 278