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 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; 42a497884dSSatish Balay PetscInt buf[4096],*irowm=NULL,*icolm; /* suppress maybe-uninitialized warning */ 432c0dbf93SJed Brown 442c0dbf93SJed Brown PetscFunctionBegin; 452c0dbf93SJed Brown if (!nrow || !ncol) PetscFunctionReturn(0); 462c0dbf93SJed Brown IndexSpaceGet(buf,nrow,ncol,irowm,icolm); 4745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr); 4845b6f7e9SBarry 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 static PetscErrorCode MatSetValuesBlockedLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 552c0dbf93SJed Brown { 562c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef*)A->data; 572c0dbf93SJed Brown PetscErrorCode ierr; 586e191704SLawrence Mitchell PetscInt rbs,cbs,buf[4096],*irowm,*icolm; 592c0dbf93SJed Brown 602c0dbf93SJed Brown PetscFunctionBegin; 616e191704SLawrence Mitchell ierr = MatGetBlockSizes(A,&rbs,&cbs);CHKERRQ(ierr); 626e191704SLawrence Mitchell IndexSpaceGet(buf,nrow*rbs,ncol*cbs,irowm,icolm); 636e191704SLawrence Mitchell BlockIndicesExpand(nrow,irow,rbs,irowm); 646e191704SLawrence Mitchell BlockIndicesExpand(ncol,icol,cbs,icolm); 656e191704SLawrence Mitchell ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow*rbs,irowm,irowm);CHKERRQ(ierr); 666e191704SLawrence Mitchell ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol*cbs,icolm,icolm);CHKERRQ(ierr); 676e191704SLawrence Mitchell ierr = (*lr->SetValues)(lr->Top,nrow*rbs,irowm,ncol*cbs,icolm,y,addv);CHKERRQ(ierr); 686e191704SLawrence Mitchell IndexSpaceRestore(buf,nrow*rbs,ncol*cbs,irowm,icolm); 692c0dbf93SJed Brown PetscFunctionReturn(0); 702c0dbf93SJed Brown } 712c0dbf93SJed Brown 722c0dbf93SJed Brown static PetscErrorCode MatSetValuesLocal_LocalRef_Scalar(Mat A,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv) 732c0dbf93SJed Brown { 742c0dbf93SJed Brown Mat_LocalRef *lr = (Mat_LocalRef*)A->data; 752c0dbf93SJed Brown PetscErrorCode ierr; 762c0dbf93SJed Brown PetscInt buf[4096],*irowm,*icolm; 772c0dbf93SJed Brown 782c0dbf93SJed Brown PetscFunctionBegin; 792c0dbf93SJed Brown IndexSpaceGet(buf,nrow,ncol,irowm,icolm); 8089f730bfSLawrence Mitchell /* If the row IS defining this submatrix was an ISBLOCK, then the unblocked LGMapApply is the right one to use. If 8189f730bfSLawrence Mitchell * instead it was (say) an ISSTRIDE with a block size > 1, then we need to use LGMapApplyBlock */ 8289f730bfSLawrence Mitchell if (lr->rowisblock) { 8389f730bfSLawrence Mitchell ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr); 8489f730bfSLawrence Mitchell } else { 8512e05225SLawrence Mitchell ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,nrow,irow,irowm);CHKERRQ(ierr); 8689f730bfSLawrence Mitchell } 8789f730bfSLawrence Mitchell /* As above, but for the column IS. */ 8889f730bfSLawrence Mitchell if (lr->colisblock) { 8989f730bfSLawrence Mitchell ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,ncol,icol,icolm);CHKERRQ(ierr); 9089f730bfSLawrence Mitchell } else { 9112e05225SLawrence Mitchell ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,ncol,icol,icolm);CHKERRQ(ierr); 9289f730bfSLawrence Mitchell } 932c0dbf93SJed Brown ierr = (*lr->SetValues)(lr->Top,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr); 942c0dbf93SJed Brown IndexSpaceRestore(buf,nrow,ncol,irowm,icolm); 952c0dbf93SJed Brown PetscFunctionReturn(0); 962c0dbf93SJed Brown } 972c0dbf93SJed Brown 982265a299SJed Brown /* Compose an IS with an ISLocalToGlobalMapping to map from IS source indices to global indices */ 992c0dbf93SJed Brown static PetscErrorCode ISL2GCompose(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog) 1002c0dbf93SJed Brown { 1012c0dbf93SJed Brown PetscErrorCode ierr; 1022c0dbf93SJed Brown const PetscInt *idx; 1032c0dbf93SJed Brown PetscInt m,*idxm; 10412e05225SLawrence Mitchell PetscInt bs; 105565245c5SBarry Smith PetscBool isblock; 1062c0dbf93SJed Brown 1072c0dbf93SJed Brown PetscFunctionBegin; 108b3e8af9fSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,1); 109b3e8af9fSJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2); 110b3e8af9fSJed Brown PetscValidPointer(cltog,3); 111565245c5SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr); 112565245c5SBarry Smith if (isblock) { 11312e05225SLawrence Mitchell PetscInt lbs; 114565245c5SBarry Smith 115565245c5SBarry Smith ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr); 116565245c5SBarry Smith ierr = ISLocalToGlobalMappingGetBlockSize(ltog,&lbs);CHKERRQ(ierr); 117565245c5SBarry Smith if (bs == lbs) { 118565245c5SBarry Smith ierr = ISGetLocalSize(is,&m);CHKERRQ(ierr); 119565245c5SBarry Smith m = m/bs; 120565245c5SBarry Smith ierr = ISBlockGetIndices(is,&idx);CHKERRQ(ierr); 121565245c5SBarry Smith ierr = PetscMalloc1(m,&idxm);CHKERRQ(ierr); 122565245c5SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm);CHKERRQ(ierr); 123565245c5SBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr); 124565245c5SBarry Smith ierr = ISBlockRestoreIndices(is,&idx);CHKERRQ(ierr); 125565245c5SBarry Smith PetscFunctionReturn(0); 126565245c5SBarry Smith } 127565245c5SBarry Smith } 1282c0dbf93SJed Brown ierr = ISGetLocalSize(is,&m);CHKERRQ(ierr); 1292c0dbf93SJed Brown ierr = ISGetIndices(is,&idx);CHKERRQ(ierr); 13012e05225SLawrence Mitchell ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr); 131785e854fSJed Brown ierr = PetscMalloc1(m,&idxm);CHKERRQ(ierr); 1322c0dbf93SJed Brown if (ltog) { 1332c0dbf93SJed Brown ierr = ISLocalToGlobalMappingApply(ltog,m,idx,idxm);CHKERRQ(ierr); 1342c0dbf93SJed Brown } else { 135580bdb30SBarry Smith ierr = PetscArraycpy(idxm,idx,m);CHKERRQ(ierr); 1362c0dbf93SJed Brown } 13712e05225SLawrence Mitchell ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr); 1382c0dbf93SJed Brown ierr = ISRestoreIndices(is,&idx);CHKERRQ(ierr); 1392c0dbf93SJed Brown PetscFunctionReturn(0); 1402c0dbf93SJed Brown } 1412c0dbf93SJed Brown 1422c0dbf93SJed Brown static PetscErrorCode ISL2GComposeBlock(IS is,ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping *cltog) 1432c0dbf93SJed Brown { 1442c0dbf93SJed Brown PetscErrorCode ierr; 1452c0dbf93SJed Brown const PetscInt *idx; 14645b6f7e9SBarry Smith PetscInt m,*idxm,bs; 1472c0dbf93SJed Brown 1482c0dbf93SJed Brown PetscFunctionBegin; 149b3e8af9fSJed Brown PetscValidHeaderSpecific(is,IS_CLASSID,1); 150b3e8af9fSJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,2); 151b3e8af9fSJed Brown PetscValidPointer(cltog,3); 1522c0dbf93SJed Brown ierr = ISBlockGetLocalSize(is,&m);CHKERRQ(ierr); 1532c0dbf93SJed Brown ierr = ISBlockGetIndices(is,&idx);CHKERRQ(ierr); 15445b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingGetBlockSize(ltog,&bs);CHKERRQ(ierr); 155785e854fSJed Brown ierr = PetscMalloc1(m,&idxm);CHKERRQ(ierr); 1562c0dbf93SJed Brown if (ltog) { 15745b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingApplyBlock(ltog,m,idx,idxm);CHKERRQ(ierr); 1582c0dbf93SJed Brown } else { 159580bdb30SBarry Smith ierr = PetscArraycpy(idxm,idx,m);CHKERRQ(ierr); 1602c0dbf93SJed Brown } 16145b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)is),bs,m,idxm,PETSC_OWN_POINTER,cltog);CHKERRQ(ierr); 1622c0dbf93SJed Brown ierr = ISBlockRestoreIndices(is,&idx);CHKERRQ(ierr); 1632c0dbf93SJed Brown PetscFunctionReturn(0); 1642c0dbf93SJed Brown } 1652c0dbf93SJed Brown 1662c0dbf93SJed Brown static PetscErrorCode MatDestroy_LocalRef(Mat B) 1672c0dbf93SJed Brown { 1682c0dbf93SJed Brown PetscErrorCode ierr; 1692c0dbf93SJed Brown 1702c0dbf93SJed Brown PetscFunctionBegin; 1712c0dbf93SJed Brown ierr = PetscFree(B->data);CHKERRQ(ierr); 1722c0dbf93SJed Brown PetscFunctionReturn(0); 1732c0dbf93SJed Brown } 1742c0dbf93SJed Brown 1752c0dbf93SJed Brown /*@ 1762c0dbf93SJed Brown MatCreateLocalRef - Gets a logical reference to a local submatrix, for use in assembly 1772c0dbf93SJed Brown 1782c0dbf93SJed Brown Not Collective 1792c0dbf93SJed Brown 1804165533cSJose E. Roman Input Parameters: 1812c0dbf93SJed Brown + A - Full matrix, generally parallel 1822c0dbf93SJed Brown . isrow - Local index set for the rows 1832c0dbf93SJed Brown - iscol - Local index set for the columns 1842c0dbf93SJed Brown 1854165533cSJose E. Roman Output Parameter: 1862c0dbf93SJed Brown . newmat - New serial Mat 1872c0dbf93SJed Brown 1882c0dbf93SJed Brown Level: developer 1892c0dbf93SJed Brown 1902c0dbf93SJed Brown Notes: 1912c0dbf93SJed Brown Most will use MatGetLocalSubMatrix() which returns a real matrix corresponding to the local 1922c0dbf93SJed Brown block if it available, such as with matrix formats that store these blocks separately. 1932c0dbf93SJed Brown 1942c0dbf93SJed Brown The new matrix forwards MatSetValuesLocal() and MatSetValuesBlockedLocal() to the global system. 1952c0dbf93SJed Brown In general, it does not define MatMult() or any other functions. Local submatrices can be nested. 1962c0dbf93SJed Brown 1972c0dbf93SJed Brown .seealso: MatSetValuesLocal(), MatSetValuesBlockedLocal(), MatGetLocalSubMatrix(), MatCreateSubMatrix() 1982c0dbf93SJed Brown @*/ 1997087cfbeSBarry Smith PetscErrorCode MatCreateLocalRef(Mat A,IS isrow,IS iscol,Mat *newmat) 2002c0dbf93SJed Brown { 2012c0dbf93SJed Brown PetscErrorCode ierr; 2022c0dbf93SJed Brown Mat_LocalRef *lr; 2032c0dbf93SJed Brown Mat B; 2042c0dbf93SJed Brown PetscInt m,n; 2052c0dbf93SJed Brown PetscBool islr; 2062c0dbf93SJed Brown 2072c0dbf93SJed Brown PetscFunctionBegin; 2082c0dbf93SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2092c0dbf93SJed Brown PetscValidHeaderSpecific(isrow,IS_CLASSID,2); 2102c0dbf93SJed Brown PetscValidHeaderSpecific(iscol,IS_CLASSID,3); 2112c0dbf93SJed Brown PetscValidPointer(newmat,4); 212*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!A->rmap->mapping,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Matrix must have local to global mapping provided before this call"); 213f4259b30SLisandro Dalcin *newmat = NULL; 2142c0dbf93SJed Brown 2152c0dbf93SJed Brown ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr); 2162c0dbf93SJed Brown ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr); 2172c0dbf93SJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 2182c0dbf93SJed Brown ierr = MatSetSizes(B,m,n,m,n);CHKERRQ(ierr); 2192c0dbf93SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)B,MATLOCALREF);CHKERRQ(ierr); 220be7c243fSJed Brown ierr = MatSetUp(B);CHKERRQ(ierr); 2212c0dbf93SJed Brown 2222c0dbf93SJed Brown B->ops->destroy = MatDestroy_LocalRef; 2232c0dbf93SJed Brown 224b00a9115SJed Brown ierr = PetscNewLog(B,&lr);CHKERRQ(ierr); 2252c0dbf93SJed Brown B->data = (void*)lr; 2262c0dbf93SJed Brown 227251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATLOCALREF,&islr);CHKERRQ(ierr); 2282265a299SJed Brown if (islr) { 2292265a299SJed Brown Mat_LocalRef *alr = (Mat_LocalRef*)A->data; 2302265a299SJed Brown lr->Top = alr->Top; 2312265a299SJed Brown } else { 2322265a299SJed Brown /* This does not increase the reference count because MatLocalRef is not allowed to live longer than its parent */ 2332265a299SJed Brown lr->Top = A; 2342265a299SJed Brown } 2352265a299SJed Brown { 2362c0dbf93SJed Brown ISLocalToGlobalMapping rltog,cltog; 2376e191704SLawrence Mitchell PetscInt arbs,acbs,rbs,cbs; 2382c0dbf93SJed Brown 2392265a299SJed Brown /* We will translate directly to global indices for the top level */ 2402c0dbf93SJed Brown lr->SetValues = MatSetValues; 2412c0dbf93SJed Brown lr->SetValuesBlocked = MatSetValuesBlocked; 2422c0dbf93SJed Brown 2432c0dbf93SJed Brown B->ops->setvalueslocal = MatSetValuesLocal_LocalRef_Scalar; 2442205254eSKarl Rupp 245992144d0SBarry Smith ierr = ISL2GCompose(isrow,A->rmap->mapping,&rltog);CHKERRQ(ierr); 246992144d0SBarry Smith if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) { 2472c0dbf93SJed Brown ierr = PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr); 2482c0dbf93SJed Brown cltog = rltog; 2492c0dbf93SJed Brown } else { 250992144d0SBarry Smith ierr = ISL2GCompose(iscol,A->cmap->mapping,&cltog);CHKERRQ(ierr); 2512c0dbf93SJed Brown } 25289f730bfSLawrence Mitchell /* Remember if the ISes we used to pull out the submatrix are of type ISBLOCK. Will be used later in 25389f730bfSLawrence Mitchell * MatSetValuesLocal_LocalRef_Scalar. */ 25489f730bfSLawrence Mitchell ierr = PetscObjectTypeCompare((PetscObject)isrow,ISBLOCK,&lr->rowisblock);CHKERRQ(ierr); 25589f730bfSLawrence Mitchell ierr = PetscObjectTypeCompare((PetscObject)iscol,ISBLOCK,&lr->colisblock);CHKERRQ(ierr); 2562c0dbf93SJed Brown ierr = MatSetLocalToGlobalMapping(B,rltog,cltog);CHKERRQ(ierr); 2576bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&rltog);CHKERRQ(ierr); 2586bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&cltog);CHKERRQ(ierr); 2592c0dbf93SJed Brown 2606e191704SLawrence Mitchell ierr = MatGetBlockSizes(A,&arbs,&acbs);CHKERRQ(ierr); 261520db06cSJed Brown ierr = ISGetBlockSize(isrow,&rbs);CHKERRQ(ierr); 262520db06cSJed Brown ierr = ISGetBlockSize(iscol,&cbs);CHKERRQ(ierr); 2636e191704SLawrence Mitchell /* Always support block interface insertion on submatrix */ 26492ea8381SJed Brown ierr = PetscLayoutSetBlockSize(B->rmap,rbs);CHKERRQ(ierr); 26592ea8381SJed Brown ierr = PetscLayoutSetBlockSize(B->cmap,cbs);CHKERRQ(ierr); 2666e191704SLawrence Mitchell if (arbs != rbs || acbs != cbs || (arbs == 1 && acbs == 1)) { 2672c0dbf93SJed Brown /* Top-level matrix has different block size, so we have to call its scalar insertion interface */ 2682c0dbf93SJed Brown B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Scalar; 2692c0dbf93SJed Brown } else { 2702c0dbf93SJed Brown /* Block sizes match so we can forward values to the top level using the block interface */ 2712c0dbf93SJed Brown B->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_LocalRef_Block; 27226fbe8dcSKarl Rupp 27345b6f7e9SBarry Smith ierr = ISL2GComposeBlock(isrow,A->rmap->mapping,&rltog);CHKERRQ(ierr); 27445b6f7e9SBarry Smith if (isrow == iscol && A->rmap->mapping == A->cmap->mapping) { 2752c0dbf93SJed Brown ierr = PetscObjectReference((PetscObject)rltog);CHKERRQ(ierr); 2762c0dbf93SJed Brown cltog = rltog; 2772c0dbf93SJed Brown } else { 27845b6f7e9SBarry Smith ierr = ISL2GComposeBlock(iscol,A->cmap->mapping,&cltog);CHKERRQ(ierr); 2792c0dbf93SJed Brown } 28045b6f7e9SBarry Smith ierr = MatSetLocalToGlobalMapping(B,rltog,cltog);CHKERRQ(ierr); 2816bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&rltog);CHKERRQ(ierr); 2826bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&cltog);CHKERRQ(ierr); 2832c0dbf93SJed Brown } 2842c0dbf93SJed Brown } 2852c0dbf93SJed Brown *newmat = B; 2862c0dbf93SJed Brown PetscFunctionReturn(0); 2872c0dbf93SJed Brown } 288