xref: /petsc/src/mat/impls/localref/mlocalref.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 { \
25*48a46eb9SPierre Jolivet     if (nrow + ncol > (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) PetscCall(PetscFree2(irowm, icolm)); \
262c0dbf93SJed Brown   } while (0)
272c0dbf93SJed Brown 
289371c9d4SSatish Balay static void BlockIndicesExpand(PetscInt n, const PetscInt idx[], PetscInt bs, PetscInt idxm[]) {
292c0dbf93SJed Brown   PetscInt i, j;
302c0dbf93SJed Brown   for (i = 0; i < n; i++) {
319371c9d4SSatish Balay     for (j = 0; j < bs; j++) { idxm[i * bs + j] = idx[i] * bs + j; }
322c0dbf93SJed Brown   }
332c0dbf93SJed Brown }
342c0dbf93SJed Brown 
359371c9d4SSatish Balay static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Block(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) {
362c0dbf93SJed Brown   Mat_LocalRef *lr = (Mat_LocalRef *)A->data;
37a497884dSSatish Balay   PetscInt      buf[4096], *irowm = NULL, *icolm; /* suppress maybe-uninitialized warning */
382c0dbf93SJed Brown 
392c0dbf93SJed Brown   PetscFunctionBegin;
402c0dbf93SJed Brown   if (!nrow || !ncol) PetscFunctionReturn(0);
412c0dbf93SJed Brown   IndexSpaceGet(buf, nrow, ncol, irowm, icolm);
429566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow, irow, irowm));
439566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol, icol, icolm));
449566063dSJacob Faibussowitsch   PetscCall((*lr->SetValuesBlocked)(lr->Top, nrow, irowm, ncol, icolm, y, addv));
452c0dbf93SJed Brown   IndexSpaceRestore(buf, nrow, ncol, irowm, icolm);
462c0dbf93SJed Brown   PetscFunctionReturn(0);
472c0dbf93SJed Brown }
482c0dbf93SJed Brown 
499371c9d4SSatish Balay static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Scalar(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) {
502c0dbf93SJed Brown   Mat_LocalRef *lr = (Mat_LocalRef *)A->data;
516e191704SLawrence Mitchell   PetscInt      rbs, cbs, buf[4096], *irowm, *icolm;
522c0dbf93SJed Brown 
532c0dbf93SJed Brown   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A, &rbs, &cbs));
556e191704SLawrence Mitchell   IndexSpaceGet(buf, nrow * rbs, ncol * cbs, irowm, icolm);
566e191704SLawrence Mitchell   BlockIndicesExpand(nrow, irow, rbs, irowm);
576e191704SLawrence Mitchell   BlockIndicesExpand(ncol, icol, cbs, icolm);
589566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow * rbs, irowm, irowm));
599566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol * cbs, icolm, icolm));
609566063dSJacob Faibussowitsch   PetscCall((*lr->SetValues)(lr->Top, nrow * rbs, irowm, ncol * cbs, icolm, y, addv));
616e191704SLawrence Mitchell   IndexSpaceRestore(buf, nrow * rbs, ncol * cbs, irowm, icolm);
622c0dbf93SJed Brown   PetscFunctionReturn(0);
632c0dbf93SJed Brown }
642c0dbf93SJed Brown 
659371c9d4SSatish Balay static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) {
662c0dbf93SJed Brown   Mat_LocalRef *lr = (Mat_LocalRef *)A->data;
672c0dbf93SJed Brown   PetscInt      buf[4096], *irowm, *icolm;
682c0dbf93SJed Brown 
692c0dbf93SJed Brown   PetscFunctionBegin;
702c0dbf93SJed Brown   IndexSpaceGet(buf, nrow, ncol, irowm, icolm);
7189f730bfSLawrence Mitchell   /* If the row IS defining this submatrix was an ISBLOCK, then the unblocked LGMapApply is the right one to use.  If
7289f730bfSLawrence Mitchell    * instead it was (say) an ISSTRIDE with a block size > 1, then we need to use LGMapApplyBlock */
7389f730bfSLawrence Mitchell   if (lr->rowisblock) {
749566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(A->rmap->mapping, nrow, irow, irowm));
7589f730bfSLawrence Mitchell   } else {
769566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping, nrow, irow, irowm));
7789f730bfSLawrence Mitchell   }
7889f730bfSLawrence Mitchell   /* As above, but for the column IS. */
7989f730bfSLawrence Mitchell   if (lr->colisblock) {
809566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(A->cmap->mapping, ncol, icol, icolm));
8189f730bfSLawrence Mitchell   } else {
829566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping, ncol, icol, icolm));
8389f730bfSLawrence Mitchell   }
849566063dSJacob Faibussowitsch   PetscCall((*lr->SetValues)(lr->Top, nrow, irowm, ncol, icolm, y, addv));
852c0dbf93SJed Brown   IndexSpaceRestore(buf, nrow, ncol, irowm, icolm);
862c0dbf93SJed Brown   PetscFunctionReturn(0);
872c0dbf93SJed Brown }
882c0dbf93SJed Brown 
892265a299SJed Brown /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */
909371c9d4SSatish Balay static PetscErrorCode ISL2GCompose(IS is, ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *cltog) {
912c0dbf93SJed Brown   const PetscInt *idx;
922c0dbf93SJed Brown   PetscInt        m, *idxm;
9312e05225SLawrence Mitchell   PetscInt        bs;
94565245c5SBarry Smith   PetscBool       isblock;
952c0dbf93SJed Brown 
962c0dbf93SJed Brown   PetscFunctionBegin;
97b3e8af9fSJed Brown   PetscValidHeaderSpecific(is, IS_CLASSID, 1);
98b3e8af9fSJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 2);
99b3e8af9fSJed Brown   PetscValidPointer(cltog, 3);
1009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)is, ISBLOCK, &isblock));
101565245c5SBarry Smith   if (isblock) {
10212e05225SLawrence Mitchell     PetscInt lbs;
103565245c5SBarry Smith 
1049566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(is, &bs));
1059566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog, &lbs));
106565245c5SBarry Smith     if (bs == lbs) {
1079566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(is, &m));
108565245c5SBarry Smith       m = m / bs;
1099566063dSJacob Faibussowitsch       PetscCall(ISBlockGetIndices(is, &idx));
1109566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m, &idxm));
1119566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingApplyBlock(ltog, m, idx, idxm));
1129566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog));
1139566063dSJacob Faibussowitsch       PetscCall(ISBlockRestoreIndices(is, &idx));
114565245c5SBarry Smith       PetscFunctionReturn(0);
115565245c5SBarry Smith     }
116565245c5SBarry Smith   }
1179566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &m));
1189566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &idx));
1199566063dSJacob Faibussowitsch   PetscCall(ISGetBlockSize(is, &bs));
1209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &idxm));
1212c0dbf93SJed Brown   if (ltog) {
1229566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(ltog, m, idx, idxm));
1232c0dbf93SJed Brown   } else {
1249566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(idxm, idx, m));
1252c0dbf93SJed Brown   }
1269566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog));
1279566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &idx));
1282c0dbf93SJed Brown   PetscFunctionReturn(0);
1292c0dbf93SJed Brown }
1302c0dbf93SJed Brown 
1319371c9d4SSatish Balay static PetscErrorCode ISL2GComposeBlock(IS is, ISLocalToGlobalMapping ltog, ISLocalToGlobalMapping *cltog) {
1322c0dbf93SJed Brown   const PetscInt *idx;
13345b6f7e9SBarry Smith   PetscInt        m, *idxm, bs;
1342c0dbf93SJed Brown 
1352c0dbf93SJed Brown   PetscFunctionBegin;
136b3e8af9fSJed Brown   PetscValidHeaderSpecific(is, IS_CLASSID, 1);
137b3e8af9fSJed Brown   PetscValidHeaderSpecific(ltog, IS_LTOGM_CLASSID, 2);
138b3e8af9fSJed Brown   PetscValidPointer(cltog, 3);
1399566063dSJacob Faibussowitsch   PetscCall(ISBlockGetLocalSize(is, &m));
1409566063dSJacob Faibussowitsch   PetscCall(ISBlockGetIndices(is, &idx));
1419566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog, &bs));
1429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &idxm));
1432c0dbf93SJed Brown   if (ltog) {
1449566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyBlock(ltog, m, idx, idxm));
1452c0dbf93SJed Brown   } else {
1469566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(idxm, idx, m));
1472c0dbf93SJed Brown   }
1489566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is), bs, m, idxm, PETSC_OWN_POINTER, cltog));
1499566063dSJacob Faibussowitsch   PetscCall(ISBlockRestoreIndices(is, &idx));
1502c0dbf93SJed Brown   PetscFunctionReturn(0);
1512c0dbf93SJed Brown }
1522c0dbf93SJed Brown 
1539371c9d4SSatish Balay static PetscErrorCode MatDestroy_LocalRef(Mat B) {
1542c0dbf93SJed Brown   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->data));
1562c0dbf93SJed Brown   PetscFunctionReturn(0);
1572c0dbf93SJed Brown }
1582c0dbf93SJed Brown 
1592c0dbf93SJed Brown /*@
1602c0dbf93SJed Brown    MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly
1612c0dbf93SJed Brown 
1622c0dbf93SJed Brown    Not Collective
1632c0dbf93SJed Brown 
1644165533cSJose E. Roman    Input Parameters:
1652c0dbf93SJed Brown + A - Full matrix, generally parallel
1662c0dbf93SJed Brown . isrow - Local index set for the rows
1672c0dbf93SJed Brown - iscol - Local index set for the columns
1682c0dbf93SJed Brown 
1694165533cSJose E. Roman    Output Parameter:
1702c0dbf93SJed Brown . newmat - New serial Mat
1712c0dbf93SJed Brown 
1722c0dbf93SJed Brown    Level: developer
1732c0dbf93SJed Brown 
1742c0dbf93SJed Brown    Notes:
1752c0dbf93SJed Brown    Most will use MatGetLocalSubMatrix() which returns a real matrix corresponding to the local
1762c0dbf93SJed Brown    block if it available, such as with matrix formats that store these blocks separately.
1772c0dbf93SJed Brown 
1782c0dbf93SJed Brown    The new matrix forwards MatSetValuesLocal() and MatSetValuesBlockedLocal() to the global system.
1792c0dbf93SJed Brown    In general, it does not define MatMult() or any other functions.  Local submatrices can be nested.
1802c0dbf93SJed Brown 
181db781477SPatrick Sanan .seealso: `MatSetValuesLocal()`, `MatSetValuesBlockedLocal()`, `MatGetLocalSubMatrix()`, `MatCreateSubMatrix()`
1822c0dbf93SJed Brown @*/
1839371c9d4SSatish Balay PetscErrorCode MatCreateLocalRef(Mat A, IS isrow, IS iscol, Mat *newmat) {
1842c0dbf93SJed Brown   Mat_LocalRef *lr;
1852c0dbf93SJed Brown   Mat           B;
1862c0dbf93SJed Brown   PetscInt      m, n;
1872c0dbf93SJed Brown   PetscBool     islr;
1882c0dbf93SJed Brown 
1892c0dbf93SJed Brown   PetscFunctionBegin;
1902c0dbf93SJed Brown   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1912c0dbf93SJed Brown   PetscValidHeaderSpecific(isrow, IS_CLASSID, 2);
1922c0dbf93SJed Brown   PetscValidHeaderSpecific(iscol, IS_CLASSID, 3);
1932c0dbf93SJed Brown   PetscValidPointer(newmat, 4);
19428b400f6SJacob Faibussowitsch   PetscCheck(A->rmap->mapping, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Matrix must have local to global mapping provided before this call");
195f4259b30SLisandro Dalcin   *newmat = NULL;
1962c0dbf93SJed Brown 
1979566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &B));
1989566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrow, &m));
1999566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscol, &n));
2009566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, m, n, m, n));
2019566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATLOCALREF));
2029566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
2032c0dbf93SJed Brown 
2042c0dbf93SJed Brown   B->ops->destroy = MatDestroy_LocalRef;
2052c0dbf93SJed Brown 
2069566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B, &lr));
2072c0dbf93SJed Brown   B->data = (void *)lr;
2082c0dbf93SJed Brown 
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATLOCALREF, &islr));
2102265a299SJed Brown   if (islr) {
2112265a299SJed Brown     Mat_LocalRef *alr = (Mat_LocalRef *)A->data;
2122265a299SJed Brown     lr->Top           = alr->Top;
2132265a299SJed Brown   } else {
2142265a299SJed Brown     /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */
2152265a299SJed Brown     lr->Top = A;
2162265a299SJed Brown   }
2172265a299SJed Brown   {
2182c0dbf93SJed Brown     ISLocalToGlobalMapping rltog, cltog;
2196e191704SLawrence Mitchell     PetscInt               arbs, acbs, rbs, cbs;
2202c0dbf93SJed Brown 
2212265a299SJed Brown     /* We will translate directly to global indices for the top level */
2222c0dbf93SJed Brown     lr->SetValues        = MatSetValues;
2232c0dbf93SJed Brown     lr->SetValuesBlocked = MatSetValuesBlocked;
2242c0dbf93SJed Brown 
2252c0dbf93SJed Brown     B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar;
2262205254eSKarl Rupp 
2279566063dSJacob Faibussowitsch     PetscCall(ISL2GCompose(isrow, A->rmap->mapping, &rltog));
228992144d0SBarry Smith     if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
2299566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)rltog));
2302c0dbf93SJed Brown       cltog = rltog;
2312c0dbf93SJed Brown     } else {
2329566063dSJacob Faibussowitsch       PetscCall(ISL2GCompose(iscol, A->cmap->mapping, &cltog));
2332c0dbf93SJed Brown     }
23489f730bfSLawrence Mitchell     /* Remember if the ISes we used to pull out the submatrix are of type ISBLOCK.  Will be used later in
23589f730bfSLawrence Mitchell      * MatSetValuesLocal_LocalRef_Scalar. */
2369566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)isrow, ISBLOCK, &lr->rowisblock));
2379566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)iscol, ISBLOCK, &lr->colisblock));
2389566063dSJacob Faibussowitsch     PetscCall(MatSetLocalToGlobalMapping(B, rltog, cltog));
2399566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&rltog));
2409566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&cltog));
2412c0dbf93SJed Brown 
2429566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSizes(A, &arbs, &acbs));
2439566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(isrow, &rbs));
2449566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(iscol, &cbs));
2456e191704SLawrence Mitchell     /* Always support block interface insertion on submatrix */
2469566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(B->rmap, rbs));
2479566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(B->cmap, cbs));
2486e191704SLawrence Mitchell     if (arbs != rbs || acbs != cbs || (arbs == 1 && acbs == 1)) {
2492c0dbf93SJed Brown       /* Top-level matrix has different block size, so we have to call its scalar insertion interface */
2502c0dbf93SJed Brown       B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar;
2512c0dbf93SJed Brown     } else {
2522c0dbf93SJed Brown       /* Block sizes match so we can forward values to the top level using the block interface */
2532c0dbf93SJed Brown       B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block;
25426fbe8dcSKarl Rupp 
2559566063dSJacob Faibussowitsch       PetscCall(ISL2GComposeBlock(isrow, A->rmap->mapping, &rltog));
25645b6f7e9SBarry Smith       if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
2579566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)rltog));
2582c0dbf93SJed Brown         cltog = rltog;
2592c0dbf93SJed Brown       } else {
2609566063dSJacob Faibussowitsch         PetscCall(ISL2GComposeBlock(iscol, A->cmap->mapping, &cltog));
2612c0dbf93SJed Brown       }
2629566063dSJacob Faibussowitsch       PetscCall(MatSetLocalToGlobalMapping(B, rltog, cltog));
2639566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingDestroy(&rltog));
2649566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingDestroy(&cltog));
2652c0dbf93SJed Brown     }
2662c0dbf93SJed Brown   }
2672c0dbf93SJed Brown   *newmat = B;
2682c0dbf93SJed Brown   PetscFunctionReturn(0);
2692c0dbf93SJed Brown }
270