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 { \ 14a52dc253SJed Brown if (nrow + ncol > (PetscInt)(sizeof(buf)/sizeof(buf[0]))) { \ 15dcca6d9dSJed Brown ierr = PetscMalloc2(nrow,&irowm,ncol,&icolm);CHKERRQ(ierr); \ 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 { \ 23a52dc253SJed Brown if (nrow + ncol > (PetscInt)(sizeof(buf)/sizeof(buf[0]))) { \ 242c0dbf93SJed Brown ierr = PetscFree2(irowm,icolm);CHKERRQ(ierr); \ 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 #undef __FUNCT__ 392c0dbf93SJed Brown #define __FUNCT__ "MatSetValuesBlockedLocal_LocalRef_Block" 402c0dbf93SJed Brown static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Block(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 412c0dbf93SJed Brown { 422c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef*)A->data; 432c0dbf93SJed Brown PetscErrorCode ierr; 44*a497884dSSatish Balay PetscInt buf[4096],*irowm=NULL,*icolm; /* suppress maybe-uninitialized warning */ 452c0dbf93SJed Brown 462c0dbf93SJed Brown PetscFunctionBegin; 472c0dbf93SJed Brown if (!nrow || !ncol) PetscFunctionReturn(0); 482c0dbf93SJed Brown IndexSpaceGet(buf,nrow,ncol,irowm,icolm); 4945b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr); 5045b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol,icol,icolm);CHKERRQ(ierr); 512c0dbf93SJed Brown ierr = (*lr->SetValuesBlocked)(lr->Top,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr); 522c0dbf93SJed Brown IndexSpaceRestore(buf,nrow,ncol,irowm,icolm); 532c0dbf93SJed Brown PetscFunctionReturn(0); 542c0dbf93SJed Brown } 552c0dbf93SJed Brown 562c0dbf93SJed Brown #undef __FUNCT__ 572c0dbf93SJed Brown #define __FUNCT__ "MatSetValuesBlockedLocal_LocalRef_Scalar" 582c0dbf93SJed Brown static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 592c0dbf93SJed Brown { 602c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef*)A->data; 612c0dbf93SJed Brown PetscErrorCode ierr; 626e191704SLawrence Mitchell PetscInt rbs,cbs,buf[4096],*irowm,*icolm; 632c0dbf93SJed Brown 642c0dbf93SJed Brown PetscFunctionBegin; 656e191704SLawrence Mitchell ierr = MatGetBlockSizes(A,&rbs,&cbs);CHKERRQ(ierr); 666e191704SLawrence Mitchell IndexSpaceGet(buf,nrow*rbs,ncol*cbs,irowm,icolm); 676e191704SLawrence Mitchell BlockIndicesExpand(nrow,irow,rbs,irowm); 686e191704SLawrence Mitchell BlockIndicesExpand(ncol,icol,cbs,icolm); 696e191704SLawrence Mitchell ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow*rbs,irowm,irowm);CHKERRQ(ierr); 706e191704SLawrence Mitchell ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol*cbs,icolm,icolm);CHKERRQ(ierr); 716e191704SLawrence Mitchell ierr = (*lr->SetValues)(lr->Top,nrow*rbs,irowm,ncol*cbs,icolm,y,addv);CHKERRQ(ierr); 726e191704SLawrence Mitchell IndexSpaceRestore(buf,nrow*rbs,ncol*cbs,irowm,icolm); 732c0dbf93SJed Brown PetscFunctionReturn(0); 742c0dbf93SJed Brown } 752c0dbf93SJed Brown 762c0dbf93SJed Brown #undef __FUNCT__ 772c0dbf93SJed Brown #define __FUNCT__ "MatSetValuesLocal_LocalRef_Scalar" 782c0dbf93SJed Brown static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 792c0dbf93SJed Brown { 802c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef*)A->data; 812c0dbf93SJed Brown PetscErrorCode ierr; 822c0dbf93SJed Brown PetscInt buf[4096],*irowm,*icolm; 832c0dbf93SJed Brown 842c0dbf93SJed Brown PetscFunctionBegin; 852c0dbf93SJed Brown IndexSpaceGet(buf,nrow,ncol,irowm,icolm); 8689f730bfSLawrence Mitchell /* If the row IS defining this submatrix was an ISBLOCK, then the unblocked LGMapApply is the right one to use. If 8789f730bfSLawrence Mitchell * instead it was (say) an ISSTRIDE with a block size > 1, then we need to use LGMapApplyBlock */ 8889f730bfSLawrence Mitchell if (lr->rowisblock) { 8989f730bfSLawrence Mitchell ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr); 9089f730bfSLawrence Mitchell } else { 9112e05225SLawrence Mitchell ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr); 9289f730bfSLawrence Mitchell } 9389f730bfSLawrence Mitchell /* As above, but for the column IS. */ 9489f730bfSLawrence Mitchell if (lr->colisblock) { 9589f730bfSLawrence Mitchell ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,ncol,icol,icolm);CHKERRQ(ierr); 9689f730bfSLawrence Mitchell } else { 9712e05225SLawrence Mitchell ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol,icol,icolm);CHKERRQ(ierr); 9889f730bfSLawrence Mitchell } 992c0dbf93SJed Brown ierr = (*lr->SetValues)(lr->Top,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr); 1002c0dbf93SJed Brown IndexSpaceRestore(buf,nrow,ncol,irowm,icolm); 1012c0dbf93SJed Brown PetscFunctionReturn(0); 1022c0dbf93SJed Brown } 1032c0dbf93SJed Brown 1042c0dbf93SJed Brown #undef __FUNCT__ 1052c0dbf93SJed Brown #define __FUNCT__ "ISL2GCompose" 1062265a299SJed Brown /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */ 1072c0dbf93SJed Brown static PetscErrorCode ISL2GCompose(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog) 1082c0dbf93SJed Brown { 1092c0dbf93SJed Brown PetscErrorCode ierr; 1102c0dbf93SJed Brown const PetscInt *idx; 1112c0dbf93SJed Brown PetscInt m,*idxm; 11212e05225SLawrence Mitchell PetscInt bs; 113565245c5SBarry Smith PetscBool isblock; 1142c0dbf93SJed Brown 1152c0dbf93SJed Brown PetscFunctionBegin; 116b3e8af9fSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,1); 117b3e8af9fSJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2); 118b3e8af9fSJed Brown PetscValidPointer(cltog,3); 119565245c5SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr); 120565245c5SBarry Smith if (isblock) { 12112e05225SLawrence Mitchell PetscInt lbs; 122565245c5SBarry Smith 123565245c5SBarry Smith ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr); 124565245c5SBarry Smith ierr = ISLocalToGlobalMappingGetBlockSize(ltog,&lbs);CHKERRQ(ierr); 125565245c5SBarry Smith if (bs == lbs) { 126565245c5SBarry Smith ierr = ISGetLocalSize(is,&m);CHKERRQ(ierr); 127565245c5SBarry Smith m = m/bs; 128565245c5SBarry Smith ierr = ISBlockGetIndices(is,&idx);CHKERRQ(ierr); 129565245c5SBarry Smith ierr = PetscMalloc1(m,&idxm);CHKERRQ(ierr); 130565245c5SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm);CHKERRQ(ierr); 131565245c5SBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr); 132565245c5SBarry Smith ierr = ISBlockRestoreIndices(is,&idx);CHKERRQ(ierr); 133565245c5SBarry Smith PetscFunctionReturn(0); 134565245c5SBarry Smith } 135565245c5SBarry Smith } 1362c0dbf93SJed Brown ierr = ISGetLocalSize(is,&m);CHKERRQ(ierr); 1372c0dbf93SJed Brown ierr = ISGetIndices(is,&idx);CHKERRQ(ierr); 13812e05225SLawrence Mitchell ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr); 139785e854fSJed Brown ierr = PetscMalloc1(m,&idxm);CHKERRQ(ierr); 1402c0dbf93SJed Brown if (ltog) { 1412c0dbf93SJed Brown ierr = ISLocalToGlobalMappingApply(ltog,m,idx,idxm);CHKERRQ(ierr); 1422c0dbf93SJed Brown } else { 1432c0dbf93SJed Brown ierr = PetscMemcpy(idxm,idx,m*sizeof(PetscInt));CHKERRQ(ierr); 1442c0dbf93SJed Brown } 14512e05225SLawrence Mitchell ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr); 1462c0dbf93SJed Brown ierr = ISRestoreIndices(is,&idx);CHKERRQ(ierr); 1472c0dbf93SJed Brown PetscFunctionReturn(0); 1482c0dbf93SJed Brown } 1492c0dbf93SJed Brown 1502c0dbf93SJed Brown #undef __FUNCT__ 1512c0dbf93SJed Brown #define __FUNCT__ "ISL2GComposeBlock" 1522c0dbf93SJed Brown static PetscErrorCode ISL2GComposeBlock(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog) 1532c0dbf93SJed Brown { 1542c0dbf93SJed Brown PetscErrorCode ierr; 1552c0dbf93SJed Brown const PetscInt *idx; 15645b6f7e9SBarry Smith PetscInt m,*idxm,bs; 1572c0dbf93SJed Brown 1582c0dbf93SJed Brown PetscFunctionBegin; 159b3e8af9fSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,1); 160b3e8af9fSJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2); 161b3e8af9fSJed Brown PetscValidPointer(cltog,3); 1622c0dbf93SJed Brown ierr = ISBlockGetLocalSize(is,&m);CHKERRQ(ierr); 1632c0dbf93SJed Brown ierr = ISBlockGetIndices(is,&idx);CHKERRQ(ierr); 16445b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingGetBlockSize(ltog,&bs);CHKERRQ(ierr); 165785e854fSJed Brown ierr = PetscMalloc1(m,&idxm);CHKERRQ(ierr); 1662c0dbf93SJed Brown if (ltog) { 16745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm);CHKERRQ(ierr); 1682c0dbf93SJed Brown } else { 1692c0dbf93SJed Brown ierr = PetscMemcpy(idxm,idx,m*sizeof(PetscInt));CHKERRQ(ierr); 1702c0dbf93SJed Brown } 17145b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr); 1722c0dbf93SJed Brown ierr = ISBlockRestoreIndices(is,&idx);CHKERRQ(ierr); 1732c0dbf93SJed Brown PetscFunctionReturn(0); 1742c0dbf93SJed Brown } 1752c0dbf93SJed Brown 1762c0dbf93SJed Brown #undef __FUNCT__ 1772c0dbf93SJed Brown #define __FUNCT__ "MatDestroy_LocalRef" 1782c0dbf93SJed Brown static PetscErrorCode MatDestroy_LocalRef(Mat B) 1792c0dbf93SJed Brown { 1802c0dbf93SJed Brown PetscErrorCode ierr; 1812c0dbf93SJed Brown 1822c0dbf93SJed Brown PetscFunctionBegin; 1832c0dbf93SJed Brown ierr = PetscFree(B->data);CHKERRQ(ierr); 1842c0dbf93SJed Brown PetscFunctionReturn(0); 1852c0dbf93SJed Brown } 1862c0dbf93SJed Brown 1872c0dbf93SJed Brown 1882c0dbf93SJed Brown #undef __FUNCT__ 1892c0dbf93SJed Brown #define __FUNCT__ "MatCreateLocalRef" 1902c0dbf93SJed Brown /*@ 1912c0dbf93SJed Brown MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly 1922c0dbf93SJed Brown 1932c0dbf93SJed Brown Not Collective 1942c0dbf93SJed Brown 1952c0dbf93SJed Brown Input Arguments: 1962c0dbf93SJed Brown + A - Full matrix, generally parallel 1972c0dbf93SJed Brown . isrow - Local index set for the rows 1982c0dbf93SJed Brown - iscol - Local index set for the columns 1992c0dbf93SJed Brown 2002c0dbf93SJed Brown Output Arguments: 2012c0dbf93SJed Brown . newmat - New serial Mat 2022c0dbf93SJed Brown 2032c0dbf93SJed Brown Level: developer 2042c0dbf93SJed Brown 2052c0dbf93SJed Brown Notes: 2062c0dbf93SJed Brown Most will use MatGetLocalSubMatrix() which returns a real matrix corresponding to the local 2072c0dbf93SJed Brown block if it available, such as with matrix formats that store these blocks separately. 2082c0dbf93SJed Brown 2092c0dbf93SJed Brown The new matrix forwards MatSetValuesLocal() and MatSetValuesBlockedLocal() to the global system. 2102c0dbf93SJed Brown In general, it does not define MatMult() or any other functions. Local submatrices can be nested. 2112c0dbf93SJed Brown 2122c0dbf93SJed Brown .seealso: MatSetValuesLocal(), MatSetValuesBlockedLocal(), MatGetLocalSubMatrix(), MatCreateSubMatrix() 2132c0dbf93SJed Brown @*/ 2147087cfbeSBarry Smith PetscErrorCode MatCreateLocalRef(Mat A,IS isrow,IS iscol,Mat *newmat) 2152c0dbf93SJed Brown { 2162c0dbf93SJed Brown PetscErrorCode ierr; 2172c0dbf93SJed Brown Mat_LocalRef *lr; 2182c0dbf93SJed Brown Mat B; 2192c0dbf93SJed Brown PetscInt m,n; 2202c0dbf93SJed Brown PetscBool islr; 2212c0dbf93SJed Brown 2222c0dbf93SJed Brown PetscFunctionBegin; 2232c0dbf93SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2242c0dbf93SJed Brown PetscValidHeaderSpecific(isrow,IS_CLASSID,2); 2252c0dbf93SJed Brown PetscValidHeaderSpecific(iscol,IS_CLASSID,3); 2262c0dbf93SJed Brown PetscValidPointer(newmat,4); 227de7e0af5SBarry Smith if (!A->rmap->mapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Matrix must have local to global mapping provided before this call"); 2282c0dbf93SJed Brown *newmat = 0; 2292c0dbf93SJed Brown 2302c0dbf93SJed Brown ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr); 2312c0dbf93SJed Brown ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr); 2322c0dbf93SJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 2332c0dbf93SJed Brown ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 2342c0dbf93SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)B,MATLOCALREF);CHKERRQ(ierr); 235be7c243fSJed Brown ierr = MatSetUp(B);CHKERRQ(ierr); 2362c0dbf93SJed Brown 2372c0dbf93SJed Brown B->ops->destroy = MatDestroy_LocalRef; 2382c0dbf93SJed Brown 239b00a9115SJed Brown ierr = PetscNewLog(B,&lr);CHKERRQ(ierr); 2402c0dbf93SJed Brown B->data = (void*)lr; 2412c0dbf93SJed Brown 242251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATLOCALREF,&islr);CHKERRQ(ierr); 2432265a299SJed Brown if (islr) { 2442265a299SJed Brown Mat_LocalRef *alr = (Mat_LocalRef*)A->data; 2452265a299SJed Brown lr->Top = alr->Top; 2462265a299SJed Brown } else { 2472265a299SJed Brown /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */ 2482265a299SJed Brown lr->Top = A; 2492265a299SJed Brown } 2502265a299SJed Brown { 2512c0dbf93SJed Brown ISLocalToGlobalMapping rltog,cltog; 2526e191704SLawrence Mitchell PetscInt arbs,acbs,rbs,cbs; 2532c0dbf93SJed Brown 2542265a299SJed Brown /* We will translate directly to global indices for the top level */ 2552c0dbf93SJed Brown lr->SetValues = MatSetValues; 2562c0dbf93SJed Brown lr->SetValuesBlocked = MatSetValuesBlocked; 2572c0dbf93SJed Brown 2582c0dbf93SJed Brown B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar; 2592205254eSKarl Rupp 260992144d0SBarry Smith ierr = ISL2GCompose(isrow,A->rmap->mapping,&rltog);CHKERRQ(ierr); 261992144d0SBarry Smith if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) { 2622c0dbf93SJed Brown ierr = PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr); 2632c0dbf93SJed Brown cltog = rltog; 2642c0dbf93SJed Brown } else { 265992144d0SBarry Smith ierr = ISL2GCompose(iscol,A->cmap->mapping,&cltog);CHKERRQ(ierr); 2662c0dbf93SJed Brown } 26789f730bfSLawrence Mitchell /* Remember if the ISes we used to pull out the submatrix are of type ISBLOCK. Will be used later in 26889f730bfSLawrence Mitchell * MatSetValuesLocal_LocalRef_Scalar. */ 26989f730bfSLawrence Mitchell ierr = PetscObjectTypeCompare((PetscObject)isrow,ISBLOCK,&lr->rowisblock);CHKERRQ(ierr); 27089f730bfSLawrence Mitchell ierr = PetscObjectTypeCompare((PetscObject)iscol,ISBLOCK,&lr->colisblock);CHKERRQ(ierr); 2712c0dbf93SJed Brown ierr = MatSetLocalToGlobalMapping(B,rltog,cltog);CHKERRQ(ierr); 2726bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&rltog);CHKERRQ(ierr); 2736bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&cltog);CHKERRQ(ierr); 2742c0dbf93SJed Brown 2756e191704SLawrence Mitchell ierr = MatGetBlockSizes(A,&arbs,&acbs);CHKERRQ(ierr); 276520db06cSJed Brown ierr = ISGetBlockSize(isrow,&rbs);CHKERRQ(ierr); 277520db06cSJed Brown ierr = ISGetBlockSize(iscol,&cbs);CHKERRQ(ierr); 2786e191704SLawrence Mitchell /* Always support block interface insertion on submatrix */ 27992ea8381SJed Brown ierr = PetscLayoutSetBlockSize(B->rmap,rbs);CHKERRQ(ierr); 28092ea8381SJed Brown ierr = PetscLayoutSetBlockSize(B->cmap,cbs);CHKERRQ(ierr); 2816e191704SLawrence Mitchell if (arbs != rbs || acbs != cbs || (arbs == 1 && acbs == 1)) { 2822c0dbf93SJed Brown /* Top-level matrix has different block size, so we have to call its scalar insertion interface */ 2832c0dbf93SJed Brown B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar; 2842c0dbf93SJed Brown } else { 2852c0dbf93SJed Brown /* Block sizes match so we can forward values to the top level using the block interface */ 2862c0dbf93SJed Brown B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block; 28726fbe8dcSKarl Rupp 28845b6f7e9SBarry Smith ierr = ISL2GComposeBlock(isrow,A->rmap->mapping,&rltog);CHKERRQ(ierr); 28945b6f7e9SBarry Smith if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) { 2902c0dbf93SJed Brown ierr = PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr); 2912c0dbf93SJed Brown cltog = rltog; 2922c0dbf93SJed Brown } else { 29345b6f7e9SBarry Smith ierr = ISL2GComposeBlock(iscol,A->cmap->mapping,&cltog);CHKERRQ(ierr); 2942c0dbf93SJed Brown } 29545b6f7e9SBarry Smith ierr = MatSetLocalToGlobalMapping(B,rltog,cltog);CHKERRQ(ierr); 2966bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&rltog);CHKERRQ(ierr); 2976bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&cltog);CHKERRQ(ierr); 2982c0dbf93SJed Brown } 2992c0dbf93SJed Brown } 3002c0dbf93SJed Brown *newmat = B; 3012c0dbf93SJed Brown PetscFunctionReturn(0); 3022c0dbf93SJed Brown } 303