xref: /petsc/src/mat/impls/localref/mlocalref.c (revision dcca6d9d80ebd869fe6029bd05a3aa9faafef49e)
12c0dbf93SJed Brown 
2b45d2f2cSJed Brown #include <petsc-private/matimpl.h>          /*I "petscmat.h" I*/
32c0dbf93SJed Brown 
42c0dbf93SJed Brown typedef struct {
52c0dbf93SJed Brown   Mat Top;
62c0dbf93SJed Brown   PetscErrorCode (*SetValues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
72c0dbf93SJed Brown   PetscErrorCode (*SetValuesBlocked)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
82c0dbf93SJed Brown } Mat_LocalRef;
92c0dbf93SJed Brown 
102c0dbf93SJed Brown /* These need to be macros because they use sizeof */
112c0dbf93SJed Brown #define IndexSpaceGet(buf,nrow,ncol,irowm,icolm) do {                   \
12a52dc253SJed Brown     if (nrow + ncol > (PetscInt)(sizeof(buf)/sizeof(buf[0]))) {         \
13*dcca6d9dSJed Brown       ierr = PetscMalloc2(nrow,&irowm,ncol,&icolm);CHKERRQ(ierr); \
142c0dbf93SJed Brown     } else {                                                            \
152c0dbf93SJed Brown       irowm = &buf[0];                                                  \
162c0dbf93SJed Brown       icolm = &buf[nrow];                                               \
172c0dbf93SJed Brown     }                                                                   \
182c0dbf93SJed Brown } while (0)
192c0dbf93SJed Brown 
202c0dbf93SJed Brown #define IndexSpaceRestore(buf,nrow,ncol,irowm,icolm) do {       \
21a52dc253SJed Brown     if (nrow + ncol > (PetscInt)(sizeof(buf)/sizeof(buf[0]))) { \
222c0dbf93SJed Brown       ierr = PetscFree2(irowm,icolm);CHKERRQ(ierr);             \
232c0dbf93SJed Brown     }                                                           \
242c0dbf93SJed Brown } while (0)
252c0dbf93SJed Brown 
262c0dbf93SJed Brown static void BlockIndicesExpand(PetscInt n,const PetscInt idx[],PetscInt bs,PetscInt idxm[])
272c0dbf93SJed Brown {
282c0dbf93SJed Brown   PetscInt i,j;
292c0dbf93SJed Brown   for (i=0; i<n; i++) {
302c0dbf93SJed Brown     for (j=0; j<bs; j++) {
312c0dbf93SJed Brown       idxm[i*bs+j] = idx[i]*bs + j;
322c0dbf93SJed Brown     }
332c0dbf93SJed Brown   }
342c0dbf93SJed Brown }
352c0dbf93SJed Brown 
362c0dbf93SJed Brown #undef __FUNCT__
372c0dbf93SJed Brown #define __FUNCT__ "MatSetValuesBlockedLocal_LocalRef_Block"
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;
412c0dbf93SJed Brown   PetscErrorCode ierr;
422c0dbf93SJed Brown   PetscInt       buf[4096],*irowm,*icolm;
432c0dbf93SJed Brown 
442c0dbf93SJed Brown   PetscFunctionBegin;
452c0dbf93SJed Brown   if (!nrow || !ncol) PetscFunctionReturn(0);
462c0dbf93SJed Brown   IndexSpaceGet(buf,nrow,ncol,irowm,icolm);
47992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->rmap->bmapping,nrow,irow,irowm);CHKERRQ(ierr);
48992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->cmap->bmapping,ncol,icol,icolm);CHKERRQ(ierr);
492c0dbf93SJed Brown   ierr = (*lr->SetValuesBlocked)(lr->Top,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr);
502c0dbf93SJed Brown   IndexSpaceRestore(buf,nrow,ncol,irowm,icolm);
512c0dbf93SJed Brown   PetscFunctionReturn(0);
522c0dbf93SJed Brown }
532c0dbf93SJed Brown 
542c0dbf93SJed Brown #undef __FUNCT__
552c0dbf93SJed Brown #define __FUNCT__ "MatSetValuesBlockedLocal_LocalRef_Scalar"
562c0dbf93SJed Brown static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
572c0dbf93SJed Brown {
582c0dbf93SJed Brown   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
592c0dbf93SJed Brown   PetscErrorCode ierr;
602c0dbf93SJed Brown   PetscInt       bs = A->rmap->bs,buf[4096],*irowm,*icolm;
612c0dbf93SJed Brown 
622c0dbf93SJed Brown   PetscFunctionBegin;
632c0dbf93SJed Brown   IndexSpaceGet(buf,nrow*bs,ncol*bs,irowm,icolm);
642c0dbf93SJed Brown   BlockIndicesExpand(nrow,irow,bs,irowm);
652c0dbf93SJed Brown   BlockIndicesExpand(ncol,icol,bs,icolm);
66992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,nrow*bs,irowm,irowm);CHKERRQ(ierr);
67992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,ncol*bs,icolm,icolm);CHKERRQ(ierr);
682c0dbf93SJed Brown   ierr = (*lr->SetValues)(lr->Top,nrow*bs,irowm,ncol*bs,icolm,y,addv);CHKERRQ(ierr);
692c0dbf93SJed Brown   IndexSpaceRestore(buf,nrow*bs,ncol*bs,irowm,icolm);
702c0dbf93SJed Brown   PetscFunctionReturn(0);
712c0dbf93SJed Brown }
722c0dbf93SJed Brown 
732c0dbf93SJed Brown #undef __FUNCT__
742c0dbf93SJed Brown #define __FUNCT__ "MatSetValuesLocal_LocalRef_Scalar"
752c0dbf93SJed Brown static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
762c0dbf93SJed Brown {
772c0dbf93SJed Brown   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
782c0dbf93SJed Brown   PetscErrorCode ierr;
792c0dbf93SJed Brown   PetscInt       buf[4096],*irowm,*icolm;
802c0dbf93SJed Brown 
812c0dbf93SJed Brown   PetscFunctionBegin;
822c0dbf93SJed Brown   IndexSpaceGet(buf,nrow,ncol,irowm,icolm);
83992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr);
84992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,ncol,icol,icolm);CHKERRQ(ierr);
852c0dbf93SJed Brown   ierr = (*lr->SetValues)(lr->Top,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr);
862c0dbf93SJed Brown   IndexSpaceRestore(buf,nrow,ncol,irowm,icolm);
872c0dbf93SJed Brown   PetscFunctionReturn(0);
882c0dbf93SJed Brown }
892c0dbf93SJed Brown 
902c0dbf93SJed Brown #undef __FUNCT__
912c0dbf93SJed Brown #define __FUNCT__ "ISL2GCompose"
922265a299SJed Brown /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */
932c0dbf93SJed Brown static PetscErrorCode ISL2GCompose(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
942c0dbf93SJed Brown {
952c0dbf93SJed Brown   PetscErrorCode ierr;
962c0dbf93SJed Brown   const PetscInt *idx;
972c0dbf93SJed Brown   PetscInt       m,*idxm;
982c0dbf93SJed Brown 
992c0dbf93SJed Brown   PetscFunctionBegin;
100b3e8af9fSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,1);
101b3e8af9fSJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2);
102b3e8af9fSJed Brown   PetscValidPointer(cltog,3);
1032c0dbf93SJed Brown   ierr = ISGetLocalSize(is,&m);CHKERRQ(ierr);
1042c0dbf93SJed Brown   ierr = ISGetIndices(is,&idx);CHKERRQ(ierr);
1052c0dbf93SJed Brown   ierr = PetscMalloc(m*sizeof(PetscInt),&idxm);CHKERRQ(ierr);
1062c0dbf93SJed Brown   if (ltog) {
1072c0dbf93SJed Brown     ierr = ISLocalToGlobalMappingApply(ltog,m,idx,idxm);CHKERRQ(ierr);
1082c0dbf93SJed Brown   } else {
1092c0dbf93SJed Brown     ierr = PetscMemcpy(idxm,idx,m*sizeof(PetscInt));CHKERRQ(ierr);
1102c0dbf93SJed Brown   }
111ce94432eSBarry Smith   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr);
1122c0dbf93SJed Brown   ierr = ISRestoreIndices(is,&idx);CHKERRQ(ierr);
1132c0dbf93SJed Brown   PetscFunctionReturn(0);
1142c0dbf93SJed Brown }
1152c0dbf93SJed Brown 
1162c0dbf93SJed Brown #undef __FUNCT__
1172c0dbf93SJed Brown #define __FUNCT__ "ISL2GComposeBlock"
1182c0dbf93SJed Brown static PetscErrorCode ISL2GComposeBlock(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
1192c0dbf93SJed Brown {
1202c0dbf93SJed Brown   PetscErrorCode ierr;
1212c0dbf93SJed Brown   const PetscInt *idx;
1222c0dbf93SJed Brown   PetscInt       m,*idxm;
1232c0dbf93SJed Brown 
1242c0dbf93SJed Brown   PetscFunctionBegin;
125b3e8af9fSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,1);
126b3e8af9fSJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2);
127b3e8af9fSJed Brown   PetscValidPointer(cltog,3);
1282c0dbf93SJed Brown   ierr = ISBlockGetLocalSize(is,&m);CHKERRQ(ierr);
1292c0dbf93SJed Brown   ierr = ISBlockGetIndices(is,&idx);CHKERRQ(ierr);
1302c0dbf93SJed Brown   ierr = PetscMalloc(m*sizeof(PetscInt),&idxm);CHKERRQ(ierr);
1312c0dbf93SJed Brown   if (ltog) {
1322c0dbf93SJed Brown     ierr = ISLocalToGlobalMappingApply(ltog,m,idx,idxm);CHKERRQ(ierr);
1332c0dbf93SJed Brown   } else {
1342c0dbf93SJed Brown     ierr = PetscMemcpy(idxm,idx,m*sizeof(PetscInt));CHKERRQ(ierr);
1352c0dbf93SJed Brown   }
136ce94432eSBarry Smith   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr);
1372c0dbf93SJed Brown   ierr = ISBlockRestoreIndices(is,&idx);CHKERRQ(ierr);
1382c0dbf93SJed Brown   PetscFunctionReturn(0);
1392c0dbf93SJed Brown }
1402c0dbf93SJed Brown 
1412c0dbf93SJed Brown #undef __FUNCT__
1422c0dbf93SJed Brown #define __FUNCT__ "MatDestroy_LocalRef"
1432c0dbf93SJed Brown static PetscErrorCode MatDestroy_LocalRef(Mat B)
1442c0dbf93SJed Brown {
1452c0dbf93SJed Brown   PetscErrorCode ierr;
1462c0dbf93SJed Brown 
1472c0dbf93SJed Brown   PetscFunctionBegin;
1482c0dbf93SJed Brown   ierr = PetscFree(B->data);CHKERRQ(ierr);
1492c0dbf93SJed Brown   PetscFunctionReturn(0);
1502c0dbf93SJed Brown }
1512c0dbf93SJed Brown 
1522c0dbf93SJed Brown 
1532c0dbf93SJed Brown #undef __FUNCT__
1542c0dbf93SJed Brown #define __FUNCT__ "MatCreateLocalRef"
1552c0dbf93SJed Brown /*@
1562c0dbf93SJed Brown    MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly
1572c0dbf93SJed Brown 
1582c0dbf93SJed Brown    Not Collective
1592c0dbf93SJed Brown 
1602c0dbf93SJed Brown    Input Arguments:
1612c0dbf93SJed Brown + A - Full matrix, generally parallel
1622c0dbf93SJed Brown . isrow - Local index set for the rows
1632c0dbf93SJed Brown - iscol - Local index set for the columns
1642c0dbf93SJed Brown 
1652c0dbf93SJed Brown    Output Arguments:
1662c0dbf93SJed Brown . newmat - New serial Mat
1672c0dbf93SJed Brown 
1682c0dbf93SJed Brown    Level: developer
1692c0dbf93SJed Brown 
1702c0dbf93SJed Brown    Notes:
1712c0dbf93SJed Brown    Most will use MatGetLocalSubMatrix() which returns a real matrix corresponding to the local
1722c0dbf93SJed Brown    block if it available, such as with matrix formats that store these blocks separately.
1732c0dbf93SJed Brown 
1742c0dbf93SJed Brown    The new matrix forwards MatSetValuesLocal() and MatSetValuesBlockedLocal() to the global system.
1752c0dbf93SJed Brown    In general, it does not define MatMult() or any other functions.  Local submatrices can be nested.
1762c0dbf93SJed Brown 
1772c0dbf93SJed Brown .seealso: MatSetValuesLocal(), MatSetValuesBlockedLocal(), MatGetLocalSubMatrix(), MatCreateSubMatrix()
1782c0dbf93SJed Brown @*/
1797087cfbeSBarry Smith PetscErrorCode  MatCreateLocalRef(Mat A,IS isrow,IS iscol,Mat *newmat)
1802c0dbf93SJed Brown {
1812c0dbf93SJed Brown   PetscErrorCode ierr;
1822c0dbf93SJed Brown   Mat_LocalRef   *lr;
1832c0dbf93SJed Brown   Mat            B;
1842c0dbf93SJed Brown   PetscInt       m,n;
1852c0dbf93SJed Brown   PetscBool      islr;
1862c0dbf93SJed Brown 
1872c0dbf93SJed Brown   PetscFunctionBegin;
1882c0dbf93SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1892c0dbf93SJed Brown   PetscValidHeaderSpecific(isrow,IS_CLASSID,2);
1902c0dbf93SJed Brown   PetscValidHeaderSpecific(iscol,IS_CLASSID,3);
1912c0dbf93SJed Brown   PetscValidPointer(newmat,4);
1922c0dbf93SJed Brown   *newmat = 0;
1932c0dbf93SJed Brown 
1942c0dbf93SJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr);
1952c0dbf93SJed Brown   ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr);
1962c0dbf93SJed Brown   ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
1972c0dbf93SJed Brown   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
1982c0dbf93SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)B,MATLOCALREF);CHKERRQ(ierr);
199be7c243fSJed Brown   ierr = MatSetUp(B);CHKERRQ(ierr);
2002c0dbf93SJed Brown 
2012c0dbf93SJed Brown   B->ops->destroy = MatDestroy_LocalRef;
2022c0dbf93SJed Brown 
2032c0dbf93SJed Brown   ierr    = PetscNewLog(B,Mat_LocalRef,&lr);CHKERRQ(ierr);
2042c0dbf93SJed Brown   B->data = (void*)lr;
2052c0dbf93SJed Brown 
206251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATLOCALREF,&islr);CHKERRQ(ierr);
2072265a299SJed Brown   if (islr) {
2082265a299SJed Brown     Mat_LocalRef *alr = (Mat_LocalRef*)A->data;
2092265a299SJed Brown     lr->Top = alr->Top;
2102265a299SJed Brown   } else {
2112265a299SJed Brown     /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */
2122265a299SJed Brown     lr->Top = A;
2132265a299SJed Brown   }
2142265a299SJed Brown   {
2152c0dbf93SJed Brown     ISLocalToGlobalMapping rltog,cltog;
2162c0dbf93SJed Brown     PetscInt               abs,rbs,cbs;
2172c0dbf93SJed Brown 
2182265a299SJed Brown     /* We will translate directly to global indices for the top level */
2192c0dbf93SJed Brown     lr->SetValues        = MatSetValues;
2202c0dbf93SJed Brown     lr->SetValuesBlocked = MatSetValuesBlocked;
2212c0dbf93SJed Brown 
2222c0dbf93SJed Brown     B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar;
2232205254eSKarl Rupp 
224992144d0SBarry Smith     ierr = ISL2GCompose(isrow,A->rmap->mapping,&rltog);CHKERRQ(ierr);
225992144d0SBarry Smith     if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
2262c0dbf93SJed Brown       ierr  = PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr);
2272c0dbf93SJed Brown       cltog = rltog;
2282c0dbf93SJed Brown     } else {
229992144d0SBarry Smith       ierr = ISL2GCompose(iscol,A->cmap->mapping,&cltog);CHKERRQ(ierr);
2302c0dbf93SJed Brown     }
2312c0dbf93SJed Brown     ierr = MatSetLocalToGlobalMapping(B,rltog,cltog);CHKERRQ(ierr);
2326bf464f9SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(&rltog);CHKERRQ(ierr);
2336bf464f9SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(&cltog);CHKERRQ(ierr);
2342c0dbf93SJed Brown 
2352c0dbf93SJed Brown     ierr = MatGetBlockSize(A,&abs);CHKERRQ(ierr);
236520db06cSJed Brown     ierr = ISGetBlockSize(isrow,&rbs);CHKERRQ(ierr);
237520db06cSJed Brown     ierr = ISGetBlockSize(iscol,&cbs);CHKERRQ(ierr);
2382c0dbf93SJed Brown     if (rbs == cbs) {           /* submatrix has block structure, so user can insert values with blocked interface */
23992ea8381SJed Brown       ierr = PetscLayoutSetBlockSize(B->rmap,rbs);CHKERRQ(ierr);
24092ea8381SJed Brown       ierr = PetscLayoutSetBlockSize(B->cmap,cbs);CHKERRQ(ierr);
24192ea8381SJed Brown       if (abs != rbs || abs == 1) {
2422c0dbf93SJed Brown         /* Top-level matrix has different block size, so we have to call its scalar insertion interface */
2432c0dbf93SJed Brown         B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar;
2442c0dbf93SJed Brown       } else {
2452c0dbf93SJed Brown         /* Block sizes match so we can forward values to the top level using the block interface */
2462c0dbf93SJed Brown         B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block;
24726fbe8dcSKarl Rupp 
248992144d0SBarry Smith         ierr = ISL2GComposeBlock(isrow,A->rmap->bmapping,&rltog);CHKERRQ(ierr);
249992144d0SBarry Smith         if (isrow == iscol && A->rmap->bmapping == A->cmap->bmapping) {
2502c0dbf93SJed Brown           ierr  =  PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr);
2512c0dbf93SJed Brown           cltog = rltog;
2522c0dbf93SJed Brown         } else {
253992144d0SBarry Smith           ierr = ISL2GComposeBlock(iscol,A->cmap->bmapping,&cltog);CHKERRQ(ierr);
2542c0dbf93SJed Brown         }
2552c0dbf93SJed Brown         ierr = MatSetLocalToGlobalMappingBlock(B,rltog,cltog);CHKERRQ(ierr);
2566bf464f9SBarry Smith         ierr = ISLocalToGlobalMappingDestroy(&rltog);CHKERRQ(ierr);
2576bf464f9SBarry Smith         ierr = ISLocalToGlobalMappingDestroy(&cltog);CHKERRQ(ierr);
2582c0dbf93SJed Brown       }
2592c0dbf93SJed Brown     }
2602c0dbf93SJed Brown   }
2612c0dbf93SJed Brown   *newmat = B;
2622c0dbf93SJed Brown   PetscFunctionReturn(0);
2632c0dbf93SJed Brown }
264