xref: /petsc/src/mat/impls/localref/mlocalref.c (revision dd39110ba674af879e763a5939b18b2b35b4c87f)
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 */
132c0dbf93SJed Brown #define IndexSpaceGet(buf,nrow,ncol,irowm,icolm) do {               \
14*dd39110bSPierre 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 
222c0dbf93SJed Brown #define IndexSpaceRestore(buf,nrow,ncol,irowm,icolm) do {           \
23*dd39110bSPierre Jolivet     if (nrow + ncol > (PetscInt)PETSC_STATIC_ARRAY_LENGTH(buf)) {   \
249566063dSJacob Faibussowitsch       PetscCall(PetscFree2(irowm,icolm));                           \
252c0dbf93SJed Brown     }                                                               \
262c0dbf93SJed Brown } while (0)
272c0dbf93SJed Brown 
282c0dbf93SJed Brown static void BlockIndicesExpand(PetscInt n,const PetscInt idx[],PetscInt bs,PetscInt idxm[])
292c0dbf93SJed Brown {
302c0dbf93SJed Brown   PetscInt i,j;
312c0dbf93SJed Brown   for (i=0; i<n; i++) {
322c0dbf93SJed Brown     for (j=0; j<bs; j++) {
332c0dbf93SJed Brown       idxm[i*bs+j] = idx[i]*bs + j;
342c0dbf93SJed Brown     }
352c0dbf93SJed Brown   }
362c0dbf93SJed Brown }
372c0dbf93SJed Brown 
382c0dbf93SJed Brown static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Block(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
392c0dbf93SJed Brown {
402c0dbf93SJed Brown   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
41a497884dSSatish Balay   PetscInt       buf[4096],*irowm=NULL,*icolm; /* suppress maybe-uninitialized warning */
422c0dbf93SJed Brown 
432c0dbf93SJed Brown   PetscFunctionBegin;
442c0dbf93SJed Brown   if (!nrow || !ncol) PetscFunctionReturn(0);
452c0dbf93SJed Brown   IndexSpaceGet(buf,nrow,ncol,irowm,icolm);
469566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm));
479566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol,icol,icolm));
489566063dSJacob Faibussowitsch   PetscCall((*lr->SetValuesBlocked)(lr->Top,nrow,irowm,ncol,icolm,y,addv));
492c0dbf93SJed Brown   IndexSpaceRestore(buf,nrow,ncol,irowm,icolm);
502c0dbf93SJed Brown   PetscFunctionReturn(0);
512c0dbf93SJed Brown }
522c0dbf93SJed Brown 
532c0dbf93SJed Brown static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
542c0dbf93SJed Brown {
552c0dbf93SJed Brown   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
566e191704SLawrence Mitchell   PetscInt       rbs,cbs,buf[4096],*irowm,*icolm;
572c0dbf93SJed Brown 
582c0dbf93SJed Brown   PetscFunctionBegin;
599566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(A,&rbs,&cbs));
606e191704SLawrence Mitchell   IndexSpaceGet(buf,nrow*rbs,ncol*cbs,irowm,icolm);
616e191704SLawrence Mitchell   BlockIndicesExpand(nrow,irow,rbs,irowm);
626e191704SLawrence Mitchell   BlockIndicesExpand(ncol,icol,cbs,icolm);
639566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow*rbs,irowm,irowm));
649566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol*cbs,icolm,icolm));
659566063dSJacob Faibussowitsch   PetscCall((*lr->SetValues)(lr->Top,nrow*rbs,irowm,ncol*cbs,icolm,y,addv));
666e191704SLawrence Mitchell   IndexSpaceRestore(buf,nrow*rbs,ncol*cbs,irowm,icolm);
672c0dbf93SJed Brown   PetscFunctionReturn(0);
682c0dbf93SJed Brown }
692c0dbf93SJed Brown 
702c0dbf93SJed Brown static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
712c0dbf93SJed Brown {
722c0dbf93SJed Brown   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
732c0dbf93SJed Brown   PetscInt       buf[4096],*irowm,*icolm;
742c0dbf93SJed Brown 
752c0dbf93SJed Brown   PetscFunctionBegin;
762c0dbf93SJed Brown   IndexSpaceGet(buf,nrow,ncol,irowm,icolm);
7789f730bfSLawrence Mitchell   /* If the row IS defining this submatrix was an ISBLOCK, then the unblocked LGMapApply is the right one to use.  If
7889f730bfSLawrence Mitchell    * instead it was (say) an ISSTRIDE with a block size > 1, then we need to use LGMapApplyBlock */
7989f730bfSLawrence Mitchell   if (lr->rowisblock) {
809566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(A->rmap->mapping,nrow,irow,irowm));
8189f730bfSLawrence Mitchell   } else {
829566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm));
8389f730bfSLawrence Mitchell   }
8489f730bfSLawrence Mitchell   /* As above, but for the column IS. */
8589f730bfSLawrence Mitchell   if (lr->colisblock) {
869566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(A->cmap->mapping,ncol,icol,icolm));
8789f730bfSLawrence Mitchell   } else {
889566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol,icol,icolm));
8989f730bfSLawrence Mitchell   }
909566063dSJacob Faibussowitsch   PetscCall((*lr->SetValues)(lr->Top,nrow,irowm,ncol,icolm,y,addv));
912c0dbf93SJed Brown   IndexSpaceRestore(buf,nrow,ncol,irowm,icolm);
922c0dbf93SJed Brown   PetscFunctionReturn(0);
932c0dbf93SJed Brown }
942c0dbf93SJed Brown 
952265a299SJed Brown /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */
962c0dbf93SJed Brown static PetscErrorCode ISL2GCompose(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
972c0dbf93SJed Brown {
982c0dbf93SJed Brown   const PetscInt *idx;
992c0dbf93SJed Brown   PetscInt       m,*idxm;
10012e05225SLawrence Mitchell   PetscInt       bs;
101565245c5SBarry Smith   PetscBool      isblock;
1022c0dbf93SJed Brown 
1032c0dbf93SJed Brown   PetscFunctionBegin;
104b3e8af9fSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,1);
105b3e8af9fSJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2);
106b3e8af9fSJed Brown   PetscValidPointer(cltog,3);
1079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock));
108565245c5SBarry Smith   if (isblock) {
10912e05225SLawrence Mitchell     PetscInt lbs;
110565245c5SBarry Smith 
1119566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(is,&bs));
1129566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog,&lbs));
113565245c5SBarry Smith     if (bs == lbs) {
1149566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(is,&m));
115565245c5SBarry Smith       m    = m/bs;
1169566063dSJacob Faibussowitsch       PetscCall(ISBlockGetIndices(is,&idx));
1179566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m,&idxm));
1189566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm));
1199566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog));
1209566063dSJacob Faibussowitsch       PetscCall(ISBlockRestoreIndices(is,&idx));
121565245c5SBarry Smith       PetscFunctionReturn(0);
122565245c5SBarry Smith     }
123565245c5SBarry Smith   }
1249566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is,&m));
1259566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is,&idx));
1269566063dSJacob Faibussowitsch   PetscCall(ISGetBlockSize(is,&bs));
1279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m,&idxm));
1282c0dbf93SJed Brown   if (ltog) {
1299566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApply(ltog,m,idx,idxm));
1302c0dbf93SJed Brown   } else {
1319566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(idxm,idx,m));
1322c0dbf93SJed Brown   }
1339566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog));
1349566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is,&idx));
1352c0dbf93SJed Brown   PetscFunctionReturn(0);
1362c0dbf93SJed Brown }
1372c0dbf93SJed Brown 
1382c0dbf93SJed Brown static PetscErrorCode ISL2GComposeBlock(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
1392c0dbf93SJed Brown {
1402c0dbf93SJed Brown   const PetscInt *idx;
14145b6f7e9SBarry Smith   PetscInt       m,*idxm,bs;
1422c0dbf93SJed Brown 
1432c0dbf93SJed Brown   PetscFunctionBegin;
144b3e8af9fSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,1);
145b3e8af9fSJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2);
146b3e8af9fSJed Brown   PetscValidPointer(cltog,3);
1479566063dSJacob Faibussowitsch   PetscCall(ISBlockGetLocalSize(is,&m));
1489566063dSJacob Faibussowitsch   PetscCall(ISBlockGetIndices(is,&idx));
1499566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingGetBlockSize(ltog,&bs));
1509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m,&idxm));
1512c0dbf93SJed Brown   if (ltog) {
1529566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm));
1532c0dbf93SJed Brown   } else {
1549566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(idxm,idx,m));
1552c0dbf93SJed Brown   }
1569566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog));
1579566063dSJacob Faibussowitsch   PetscCall(ISBlockRestoreIndices(is,&idx));
1582c0dbf93SJed Brown   PetscFunctionReturn(0);
1592c0dbf93SJed Brown }
1602c0dbf93SJed Brown 
1612c0dbf93SJed Brown static PetscErrorCode MatDestroy_LocalRef(Mat B)
1622c0dbf93SJed Brown {
1632c0dbf93SJed Brown   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->data));
1652c0dbf93SJed Brown   PetscFunctionReturn(0);
1662c0dbf93SJed Brown }
1672c0dbf93SJed Brown 
1682c0dbf93SJed Brown /*@
1692c0dbf93SJed Brown    MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly
1702c0dbf93SJed Brown 
1712c0dbf93SJed Brown    Not Collective
1722c0dbf93SJed Brown 
1734165533cSJose E. Roman    Input Parameters:
1742c0dbf93SJed Brown + A - Full matrix, generally parallel
1752c0dbf93SJed Brown . isrow - Local index set for the rows
1762c0dbf93SJed Brown - iscol - Local index set for the columns
1772c0dbf93SJed Brown 
1784165533cSJose E. Roman    Output Parameter:
1792c0dbf93SJed Brown . newmat - New serial Mat
1802c0dbf93SJed Brown 
1812c0dbf93SJed Brown    Level: developer
1822c0dbf93SJed Brown 
1832c0dbf93SJed Brown    Notes:
1842c0dbf93SJed Brown    Most will use MatGetLocalSubMatrix() which returns a real matrix corresponding to the local
1852c0dbf93SJed Brown    block if it available, such as with matrix formats that store these blocks separately.
1862c0dbf93SJed Brown 
1872c0dbf93SJed Brown    The new matrix forwards MatSetValuesLocal() and MatSetValuesBlockedLocal() to the global system.
1882c0dbf93SJed Brown    In general, it does not define MatMult() or any other functions.  Local submatrices can be nested.
1892c0dbf93SJed Brown 
1902c0dbf93SJed Brown .seealso: MatSetValuesLocal(), MatSetValuesBlockedLocal(), MatGetLocalSubMatrix(), MatCreateSubMatrix()
1912c0dbf93SJed Brown @*/
1927087cfbeSBarry Smith PetscErrorCode  MatCreateLocalRef(Mat A,IS isrow,IS iscol,Mat *newmat)
1932c0dbf93SJed Brown {
1942c0dbf93SJed Brown   Mat_LocalRef   *lr;
1952c0dbf93SJed Brown   Mat            B;
1962c0dbf93SJed Brown   PetscInt       m,n;
1972c0dbf93SJed Brown   PetscBool      islr;
1982c0dbf93SJed Brown 
1992c0dbf93SJed Brown   PetscFunctionBegin;
2002c0dbf93SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
2012c0dbf93SJed Brown   PetscValidHeaderSpecific(isrow,IS_CLASSID,2);
2022c0dbf93SJed Brown   PetscValidHeaderSpecific(iscol,IS_CLASSID,3);
2032c0dbf93SJed Brown   PetscValidPointer(newmat,4);
20428b400f6SJacob Faibussowitsch   PetscCheck(A->rmap->mapping,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Matrix must have local to global mapping provided before this call");
205f4259b30SLisandro Dalcin   *newmat = NULL;
2062c0dbf93SJed Brown 
2079566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&B));
2089566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrow,&m));
2099566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscol,&n));
2109566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,m,n,m,n));
2119566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B,MATLOCALREF));
2129566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
2132c0dbf93SJed Brown 
2142c0dbf93SJed Brown   B->ops->destroy = MatDestroy_LocalRef;
2152c0dbf93SJed Brown 
2169566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B,&lr));
2172c0dbf93SJed Brown   B->data = (void*)lr;
2182c0dbf93SJed Brown 
2199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATLOCALREF,&islr));
2202265a299SJed Brown   if (islr) {
2212265a299SJed Brown     Mat_LocalRef *alr = (Mat_LocalRef*)A->data;
2222265a299SJed Brown     lr->Top = alr->Top;
2232265a299SJed Brown   } else {
2242265a299SJed Brown     /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */
2252265a299SJed Brown     lr->Top = A;
2262265a299SJed Brown   }
2272265a299SJed Brown   {
2282c0dbf93SJed Brown     ISLocalToGlobalMapping rltog,cltog;
2296e191704SLawrence Mitchell     PetscInt               arbs,acbs,rbs,cbs;
2302c0dbf93SJed Brown 
2312265a299SJed Brown     /* We will translate directly to global indices for the top level */
2322c0dbf93SJed Brown     lr->SetValues        = MatSetValues;
2332c0dbf93SJed Brown     lr->SetValuesBlocked = MatSetValuesBlocked;
2342c0dbf93SJed Brown 
2352c0dbf93SJed Brown     B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar;
2362205254eSKarl Rupp 
2379566063dSJacob Faibussowitsch     PetscCall(ISL2GCompose(isrow,A->rmap->mapping,&rltog));
238992144d0SBarry Smith     if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
2399566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)rltog));
2402c0dbf93SJed Brown       cltog = rltog;
2412c0dbf93SJed Brown     } else {
2429566063dSJacob Faibussowitsch       PetscCall(ISL2GCompose(iscol,A->cmap->mapping,&cltog));
2432c0dbf93SJed Brown     }
24489f730bfSLawrence Mitchell     /* Remember if the ISes we used to pull out the submatrix are of type ISBLOCK.  Will be used later in
24589f730bfSLawrence Mitchell      * MatSetValuesLocal_LocalRef_Scalar. */
2469566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)isrow,ISBLOCK,&lr->rowisblock));
2479566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)iscol,ISBLOCK,&lr->colisblock));
2489566063dSJacob Faibussowitsch     PetscCall(MatSetLocalToGlobalMapping(B,rltog,cltog));
2499566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&rltog));
2509566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&cltog));
2512c0dbf93SJed Brown 
2529566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSizes(A,&arbs,&acbs));
2539566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(isrow,&rbs));
2549566063dSJacob Faibussowitsch     PetscCall(ISGetBlockSize(iscol,&cbs));
2556e191704SLawrence Mitchell     /* Always support block interface insertion on submatrix */
2569566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(B->rmap,rbs));
2579566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(B->cmap,cbs));
2586e191704SLawrence Mitchell     if (arbs != rbs || acbs != cbs || (arbs == 1 && acbs == 1)) {
2592c0dbf93SJed Brown       /* Top-level matrix has different block size, so we have to call its scalar insertion interface */
2602c0dbf93SJed Brown       B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar;
2612c0dbf93SJed Brown     } else {
2622c0dbf93SJed Brown       /* Block sizes match so we can forward values to the top level using the block interface */
2632c0dbf93SJed Brown       B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block;
26426fbe8dcSKarl Rupp 
2659566063dSJacob Faibussowitsch       PetscCall(ISL2GComposeBlock(isrow,A->rmap->mapping,&rltog));
26645b6f7e9SBarry Smith       if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
2679566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)rltog));
2682c0dbf93SJed Brown         cltog = rltog;
2692c0dbf93SJed Brown       } else {
2709566063dSJacob Faibussowitsch         PetscCall(ISL2GComposeBlock(iscol,A->cmap->mapping,&cltog));
2712c0dbf93SJed Brown       }
2729566063dSJacob Faibussowitsch       PetscCall(MatSetLocalToGlobalMapping(B,rltog,cltog));
2739566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingDestroy(&rltog));
2749566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingDestroy(&cltog));
2752c0dbf93SJed Brown     }
2762c0dbf93SJed Brown   }
2772c0dbf93SJed Brown   *newmat = B;
2782c0dbf93SJed Brown   PetscFunctionReturn(0);
2792c0dbf93SJed Brown }
280