xref: /petsc/src/mat/impls/localref/mlocalref.c (revision 45b6f7e9ecd564e8bb535880b270d4a8084ff211)
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]))) {         \
13dcca6d9dSJed 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);
47*45b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr);
48*45b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,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;
6033d57670SJed Brown   PetscInt       bs,buf[4096],*irowm,*icolm;
612c0dbf93SJed Brown 
622c0dbf93SJed Brown   PetscFunctionBegin;
6333d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
642c0dbf93SJed Brown   IndexSpaceGet(buf,nrow*bs,ncol*bs,irowm,icolm);
652c0dbf93SJed Brown   BlockIndicesExpand(nrow,irow,bs,irowm);
662c0dbf93SJed Brown   BlockIndicesExpand(ncol,icol,bs,icolm);
67992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,nrow*bs,irowm,irowm);CHKERRQ(ierr);
68992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,ncol*bs,icolm,icolm);CHKERRQ(ierr);
692c0dbf93SJed Brown   ierr = (*lr->SetValues)(lr->Top,nrow*bs,irowm,ncol*bs,icolm,y,addv);CHKERRQ(ierr);
702c0dbf93SJed Brown   IndexSpaceRestore(buf,nrow*bs,ncol*bs,irowm,icolm);
712c0dbf93SJed Brown   PetscFunctionReturn(0);
722c0dbf93SJed Brown }
732c0dbf93SJed Brown 
742c0dbf93SJed Brown #undef __FUNCT__
752c0dbf93SJed Brown #define __FUNCT__ "MatSetValuesLocal_LocalRef_Scalar"
762c0dbf93SJed Brown static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
772c0dbf93SJed Brown {
782c0dbf93SJed Brown   Mat_LocalRef   *lr = (Mat_LocalRef*)A->data;
792c0dbf93SJed Brown   PetscErrorCode ierr;
802c0dbf93SJed Brown   PetscInt       buf[4096],*irowm,*icolm;
812c0dbf93SJed Brown 
822c0dbf93SJed Brown   PetscFunctionBegin;
832c0dbf93SJed Brown   IndexSpaceGet(buf,nrow,ncol,irowm,icolm);
84992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr);
85992144d0SBarry Smith   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,ncol,icol,icolm);CHKERRQ(ierr);
862c0dbf93SJed Brown   ierr = (*lr->SetValues)(lr->Top,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr);
872c0dbf93SJed Brown   IndexSpaceRestore(buf,nrow,ncol,irowm,icolm);
882c0dbf93SJed Brown   PetscFunctionReturn(0);
892c0dbf93SJed Brown }
902c0dbf93SJed Brown 
912c0dbf93SJed Brown #undef __FUNCT__
922c0dbf93SJed Brown #define __FUNCT__ "ISL2GCompose"
932265a299SJed Brown /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */
942c0dbf93SJed Brown static PetscErrorCode ISL2GCompose(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
952c0dbf93SJed Brown {
962c0dbf93SJed Brown   PetscErrorCode ierr;
972c0dbf93SJed Brown   const PetscInt *idx;
982c0dbf93SJed Brown   PetscInt       m,*idxm;
992c0dbf93SJed Brown 
1002c0dbf93SJed Brown   PetscFunctionBegin;
101b3e8af9fSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,1);
102b3e8af9fSJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2);
103b3e8af9fSJed Brown   PetscValidPointer(cltog,3);
1042c0dbf93SJed Brown   ierr = ISGetLocalSize(is,&m);CHKERRQ(ierr);
1052c0dbf93SJed Brown   ierr = ISGetIndices(is,&idx);CHKERRQ(ierr);
106785e854fSJed Brown   ierr = PetscMalloc1(m,&idxm);CHKERRQ(ierr);
1072c0dbf93SJed Brown   if (ltog) {
1082c0dbf93SJed Brown     ierr = ISLocalToGlobalMappingApply(ltog,m,idx,idxm);CHKERRQ(ierr);
1092c0dbf93SJed Brown   } else {
1102c0dbf93SJed Brown     ierr = PetscMemcpy(idxm,idx,m*sizeof(PetscInt));CHKERRQ(ierr);
1112c0dbf93SJed Brown   }
112f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),1,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr);
1132c0dbf93SJed Brown   ierr = ISRestoreIndices(is,&idx);CHKERRQ(ierr);
1142c0dbf93SJed Brown   PetscFunctionReturn(0);
1152c0dbf93SJed Brown }
1162c0dbf93SJed Brown 
1172c0dbf93SJed Brown #undef __FUNCT__
1182c0dbf93SJed Brown #define __FUNCT__ "ISL2GComposeBlock"
1192c0dbf93SJed Brown static PetscErrorCode ISL2GComposeBlock(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog)
1202c0dbf93SJed Brown {
1212c0dbf93SJed Brown   PetscErrorCode ierr;
1222c0dbf93SJed Brown   const PetscInt *idx;
123*45b6f7e9SBarry Smith   PetscInt       m,*idxm,bs;
1242c0dbf93SJed Brown 
1252c0dbf93SJed Brown   PetscFunctionBegin;
126b3e8af9fSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,1);
127b3e8af9fSJed Brown   PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2);
128b3e8af9fSJed Brown   PetscValidPointer(cltog,3);
1292c0dbf93SJed Brown   ierr = ISBlockGetLocalSize(is,&m);CHKERRQ(ierr);
1302c0dbf93SJed Brown   ierr = ISBlockGetIndices(is,&idx);CHKERRQ(ierr);
131*45b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingGetBlockSize(ltog,&bs);CHKERRQ(ierr);
132785e854fSJed Brown   ierr = PetscMalloc1(m,&idxm);CHKERRQ(ierr);
1332c0dbf93SJed Brown   if (ltog) {
134*45b6f7e9SBarry Smith     ierr = ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm);CHKERRQ(ierr);
1352c0dbf93SJed Brown   } else {
1362c0dbf93SJed Brown     ierr = PetscMemcpy(idxm,idx,m*sizeof(PetscInt));CHKERRQ(ierr);
1372c0dbf93SJed Brown   }
138*45b6f7e9SBarry Smith   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr);
1392c0dbf93SJed Brown   ierr = ISBlockRestoreIndices(is,&idx);CHKERRQ(ierr);
1402c0dbf93SJed Brown   PetscFunctionReturn(0);
1412c0dbf93SJed Brown }
1422c0dbf93SJed Brown 
1432c0dbf93SJed Brown #undef __FUNCT__
1442c0dbf93SJed Brown #define __FUNCT__ "MatDestroy_LocalRef"
1452c0dbf93SJed Brown static PetscErrorCode MatDestroy_LocalRef(Mat B)
1462c0dbf93SJed Brown {
1472c0dbf93SJed Brown   PetscErrorCode ierr;
1482c0dbf93SJed Brown 
1492c0dbf93SJed Brown   PetscFunctionBegin;
1502c0dbf93SJed Brown   ierr = PetscFree(B->data);CHKERRQ(ierr);
1512c0dbf93SJed Brown   PetscFunctionReturn(0);
1522c0dbf93SJed Brown }
1532c0dbf93SJed Brown 
1542c0dbf93SJed Brown 
1552c0dbf93SJed Brown #undef __FUNCT__
1562c0dbf93SJed Brown #define __FUNCT__ "MatCreateLocalRef"
1572c0dbf93SJed Brown /*@
1582c0dbf93SJed Brown    MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly
1592c0dbf93SJed Brown 
1602c0dbf93SJed Brown    Not Collective
1612c0dbf93SJed Brown 
1622c0dbf93SJed Brown    Input Arguments:
1632c0dbf93SJed Brown + A - Full matrix, generally parallel
1642c0dbf93SJed Brown . isrow - Local index set for the rows
1652c0dbf93SJed Brown - iscol - Local index set for the columns
1662c0dbf93SJed Brown 
1672c0dbf93SJed Brown    Output Arguments:
1682c0dbf93SJed Brown . newmat - New serial Mat
1692c0dbf93SJed Brown 
1702c0dbf93SJed Brown    Level: developer
1712c0dbf93SJed Brown 
1722c0dbf93SJed Brown    Notes:
1732c0dbf93SJed Brown    Most will use MatGetLocalSubMatrix() which returns a real matrix corresponding to the local
1742c0dbf93SJed Brown    block if it available, such as with matrix formats that store these blocks separately.
1752c0dbf93SJed Brown 
1762c0dbf93SJed Brown    The new matrix forwards MatSetValuesLocal() and MatSetValuesBlockedLocal() to the global system.
1772c0dbf93SJed Brown    In general, it does not define MatMult() or any other functions.  Local submatrices can be nested.
1782c0dbf93SJed Brown 
1792c0dbf93SJed Brown .seealso: MatSetValuesLocal(), MatSetValuesBlockedLocal(), MatGetLocalSubMatrix(), MatCreateSubMatrix()
1802c0dbf93SJed Brown @*/
1817087cfbeSBarry Smith PetscErrorCode  MatCreateLocalRef(Mat A,IS isrow,IS iscol,Mat *newmat)
1822c0dbf93SJed Brown {
1832c0dbf93SJed Brown   PetscErrorCode ierr;
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);
1942c0dbf93SJed Brown   *newmat = 0;
1952c0dbf93SJed Brown 
1962c0dbf93SJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr);
1972c0dbf93SJed Brown   ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr);
1982c0dbf93SJed Brown   ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
1992c0dbf93SJed Brown   ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr);
2002c0dbf93SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)B,MATLOCALREF);CHKERRQ(ierr);
201be7c243fSJed Brown   ierr = MatSetUp(B);CHKERRQ(ierr);
2022c0dbf93SJed Brown 
2032c0dbf93SJed Brown   B->ops->destroy = MatDestroy_LocalRef;
2042c0dbf93SJed Brown 
205b00a9115SJed Brown   ierr    = PetscNewLog(B,&lr);CHKERRQ(ierr);
2062c0dbf93SJed Brown   B->data = (void*)lr;
2072c0dbf93SJed Brown 
208251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATLOCALREF,&islr);CHKERRQ(ierr);
2092265a299SJed Brown   if (islr) {
2102265a299SJed Brown     Mat_LocalRef *alr = (Mat_LocalRef*)A->data;
2112265a299SJed Brown     lr->Top = alr->Top;
2122265a299SJed Brown   } else {
2132265a299SJed Brown     /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */
2142265a299SJed Brown     lr->Top = A;
2152265a299SJed Brown   }
2162265a299SJed Brown   {
2172c0dbf93SJed Brown     ISLocalToGlobalMapping rltog,cltog;
2182c0dbf93SJed Brown     PetscInt               abs,rbs,cbs;
2192c0dbf93SJed Brown 
2202265a299SJed Brown     /* We will translate directly to global indices for the top level */
2212c0dbf93SJed Brown     lr->SetValues        = MatSetValues;
2222c0dbf93SJed Brown     lr->SetValuesBlocked = MatSetValuesBlocked;
2232c0dbf93SJed Brown 
2242c0dbf93SJed Brown     B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar;
2252205254eSKarl Rupp 
226992144d0SBarry Smith     ierr = ISL2GCompose(isrow,A->rmap->mapping,&rltog);CHKERRQ(ierr);
227992144d0SBarry Smith     if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
2282c0dbf93SJed Brown       ierr  = PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr);
2292c0dbf93SJed Brown       cltog = rltog;
2302c0dbf93SJed Brown     } else {
231992144d0SBarry Smith       ierr = ISL2GCompose(iscol,A->cmap->mapping,&cltog);CHKERRQ(ierr);
2322c0dbf93SJed Brown     }
2332c0dbf93SJed Brown     ierr = MatSetLocalToGlobalMapping(B,rltog,cltog);CHKERRQ(ierr);
2346bf464f9SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(&rltog);CHKERRQ(ierr);
2356bf464f9SBarry Smith     ierr = ISLocalToGlobalMappingDestroy(&cltog);CHKERRQ(ierr);
2362c0dbf93SJed Brown 
2372c0dbf93SJed Brown     ierr = MatGetBlockSize(A,&abs);CHKERRQ(ierr);
238520db06cSJed Brown     ierr = ISGetBlockSize(isrow,&rbs);CHKERRQ(ierr);
239520db06cSJed Brown     ierr = ISGetBlockSize(iscol,&cbs);CHKERRQ(ierr);
2402c0dbf93SJed Brown     if (rbs == cbs) {           /* submatrix has block structure, so user can insert values with blocked interface */
24192ea8381SJed Brown       ierr = PetscLayoutSetBlockSize(B->rmap,rbs);CHKERRQ(ierr);
24292ea8381SJed Brown       ierr = PetscLayoutSetBlockSize(B->cmap,cbs);CHKERRQ(ierr);
24392ea8381SJed Brown       if (abs != rbs || abs == 1) {
2442c0dbf93SJed Brown         /* Top-level matrix has different block size, so we have to call its scalar insertion interface */
2452c0dbf93SJed Brown         B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar;
2462c0dbf93SJed Brown       } else {
2472c0dbf93SJed Brown         /* Block sizes match so we can forward values to the top level using the block interface */
2482c0dbf93SJed Brown         B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block;
24926fbe8dcSKarl Rupp 
250*45b6f7e9SBarry Smith         ierr = ISL2GComposeBlock(isrow,A->rmap->mapping,&rltog);CHKERRQ(ierr);
251*45b6f7e9SBarry Smith         if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) {
2522c0dbf93SJed Brown           ierr  =  PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr);
2532c0dbf93SJed Brown           cltog = rltog;
2542c0dbf93SJed Brown         } else {
255*45b6f7e9SBarry Smith           ierr = ISL2GComposeBlock(iscol,A->cmap->mapping,&cltog);CHKERRQ(ierr);
2562c0dbf93SJed Brown         }
257*45b6f7e9SBarry Smith         ierr = MatSetLocalToGlobalMapping(B,rltog,cltog);CHKERRQ(ierr);
2586bf464f9SBarry Smith         ierr = ISLocalToGlobalMappingDestroy(&rltog);CHKERRQ(ierr);
2596bf464f9SBarry Smith         ierr = ISLocalToGlobalMappingDestroy(&cltog);CHKERRQ(ierr);
2602c0dbf93SJed Brown       }
2612c0dbf93SJed Brown     }
2622c0dbf93SJed Brown   }
2632c0dbf93SJed Brown   *newmat = B;
2642c0dbf93SJed Brown   PetscFunctionReturn(0);
2652c0dbf93SJed Brown }
266