1ee1cef2cSJed Brown 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3ee1cef2cSJed Brown 4ee1cef2cSJed Brown typedef struct { 5ee1cef2cSJed Brown IS isrow,iscol; /* rows and columns in submatrix, only used to check consistency */ 6ee1cef2cSJed Brown Vec lwork,rwork; /* work vectors inside the scatters */ 711c5f74dSStefano Zampini Vec lwork2,rwork2; /* work vectors inside the scatters */ 8ee1cef2cSJed Brown VecScatter lrestrict,rprolong; 9ee1cef2cSJed Brown Mat A; 1054e05e6cSHong Zhang } Mat_SubVirtual; 11ee1cef2cSJed Brown 12bddd1ffdSAlp Dener static PetscErrorCode MatScale_SubMatrix(Mat N,PetscScalar a) 13bddd1ffdSAlp Dener { 14bddd1ffdSAlp Dener Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 15bddd1ffdSAlp Dener 16bddd1ffdSAlp Dener PetscFunctionBegin; 175f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(Na->A,a)); 18bddd1ffdSAlp Dener PetscFunctionReturn(0); 19bddd1ffdSAlp Dener } 20bddd1ffdSAlp Dener 21bddd1ffdSAlp Dener static PetscErrorCode MatShift_SubMatrix(Mat N,PetscScalar a) 22bddd1ffdSAlp Dener { 23bddd1ffdSAlp Dener Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 24bddd1ffdSAlp Dener 25bddd1ffdSAlp Dener PetscFunctionBegin; 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(Na->A,a)); 27ee1cef2cSJed Brown PetscFunctionReturn(0); 28ee1cef2cSJed Brown } 29ee1cef2cSJed Brown 30ee1cef2cSJed Brown static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N,Vec left,Vec right) 31ee1cef2cSJed Brown { 3254e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 33ee1cef2cSJed Brown 34ee1cef2cSJed Brown PetscFunctionBegin; 35ee1cef2cSJed Brown if (right) { 365f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->rwork)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->rprolong,right,Na->rwork,INSERT_VALUES,SCATTER_FORWARD)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->rprolong,right,Na->rwork,INSERT_VALUES,SCATTER_FORWARD)); 39ee1cef2cSJed Brown } 4011c5f74dSStefano Zampini if (left) { 415f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->lwork)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->lrestrict,left,Na->lwork,INSERT_VALUES,SCATTER_REVERSE)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->lrestrict,left,Na->lwork,INSERT_VALUES,SCATTER_REVERSE)); 44ee1cef2cSJed Brown } 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(Na->A,left ? Na->lwork : NULL,right ? Na->rwork : NULL)); 4611c5f74dSStefano Zampini PetscFunctionReturn(0); 4711c5f74dSStefano Zampini } 4811c5f74dSStefano Zampini 4911c5f74dSStefano Zampini static PetscErrorCode MatGetDiagonal_SubMatrix(Mat N,Vec d) 5011c5f74dSStefano Zampini { 5111c5f74dSStefano Zampini Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 5211c5f74dSStefano Zampini 5311c5f74dSStefano Zampini PetscFunctionBegin; 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(Na->A,Na->rwork)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->rprolong,Na->rwork,d,INSERT_VALUES,SCATTER_REVERSE)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->rprolong,Na->rwork,d,INSERT_VALUES,SCATTER_REVERSE)); 57ee1cef2cSJed Brown PetscFunctionReturn(0); 58ee1cef2cSJed Brown } 59ee1cef2cSJed Brown 60ee1cef2cSJed Brown static PetscErrorCode MatMult_SubMatrix(Mat N,Vec x,Vec y) 61ee1cef2cSJed Brown { 6254e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 63ee1cef2cSJed Brown 64ee1cef2cSJed Brown PetscFunctionBegin; 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->rwork)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->rprolong,x,Na->rwork,INSERT_VALUES,SCATTER_FORWARD)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->rprolong,x,Na->rwork,INSERT_VALUES,SCATTER_FORWARD)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Na->A,Na->rwork,Na->lwork)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD)); 71ee1cef2cSJed Brown PetscFunctionReturn(0); 72ee1cef2cSJed Brown } 73ee1cef2cSJed Brown 74ee1cef2cSJed Brown static PetscErrorCode MatMultAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3) 75ee1cef2cSJed Brown { 7654e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 77ee1cef2cSJed Brown 78ee1cef2cSJed Brown PetscFunctionBegin; 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->rwork)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->rprolong,v1,Na->rwork,INSERT_VALUES,SCATTER_FORWARD)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->rprolong,v1,Na->rwork,INSERT_VALUES,SCATTER_FORWARD)); 8211c5f74dSStefano Zampini if (v1 == v2) { 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(Na->A,Na->rwork,Na->rwork,Na->lwork)); 8411c5f74dSStefano Zampini } else if (v2 == v3) { 855f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->lwork)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->lrestrict,v2,Na->lwork,INSERT_VALUES,SCATTER_REVERSE)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->lrestrict,v2,Na->lwork,INSERT_VALUES,SCATTER_REVERSE)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(Na->A,Na->rwork,Na->lwork,Na->lwork)); 891e7f0bbdSJed Brown } else { 9011c5f74dSStefano Zampini if (!Na->lwork2) { 915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(Na->lwork,&Na->lwork2)); 9211c5f74dSStefano Zampini } else { 935f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->lwork2)); 9411c5f74dSStefano Zampini } 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->lrestrict,v2,Na->lwork2,INSERT_VALUES,SCATTER_REVERSE)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->lrestrict,v2,Na->lwork2,INSERT_VALUES,SCATTER_REVERSE)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(Na->A,Na->rwork,Na->lwork2,Na->lwork)); 9811c5f74dSStefano Zampini } 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD)); 101ee1cef2cSJed Brown PetscFunctionReturn(0); 102ee1cef2cSJed Brown } 103ee1cef2cSJed Brown 104ee1cef2cSJed Brown static PetscErrorCode MatMultTranspose_SubMatrix(Mat N,Vec x,Vec y) 105ee1cef2cSJed Brown { 10654e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 107ee1cef2cSJed Brown 108ee1cef2cSJed Brown PetscFunctionBegin; 1095f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->lwork)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->lrestrict,x,Na->lwork,INSERT_VALUES,SCATTER_REVERSE)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->lrestrict,x,Na->lwork,INSERT_VALUES,SCATTER_REVERSE)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(Na->A,Na->lwork,Na->rwork)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE)); 115ee1cef2cSJed Brown PetscFunctionReturn(0); 116ee1cef2cSJed Brown } 117ee1cef2cSJed Brown 118ee1cef2cSJed Brown static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3) 119ee1cef2cSJed Brown { 12054e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 121ee1cef2cSJed Brown 122ee1cef2cSJed Brown PetscFunctionBegin; 1235f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->lwork)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->lrestrict,v1,Na->lwork,INSERT_VALUES,SCATTER_REVERSE)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->lrestrict,v1,Na->lwork,INSERT_VALUES,SCATTER_REVERSE)); 12611c5f74dSStefano Zampini if (v1 == v2) { 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(Na->A,Na->lwork,Na->lwork,Na->rwork)); 12811c5f74dSStefano Zampini } else if (v2 == v3) { 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->rwork)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->rprolong,v2,Na->rwork,INSERT_VALUES,SCATTER_FORWARD)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->rprolong,v2,Na->rwork,INSERT_VALUES,SCATTER_FORWARD)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(Na->A,Na->lwork,Na->rwork,Na->rwork)); 1331e7f0bbdSJed Brown } else { 13411c5f74dSStefano Zampini if (!Na->rwork2) { 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(Na->rwork,&Na->rwork2)); 13611c5f74dSStefano Zampini } else { 1375f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(Na->rwork2)); 13811c5f74dSStefano Zampini } 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->rprolong,v2,Na->rwork2,INSERT_VALUES,SCATTER_FORWARD)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->rprolong,v2,Na->rwork2,INSERT_VALUES,SCATTER_FORWARD)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(Na->A,Na->lwork,Na->rwork2,Na->rwork)); 14211c5f74dSStefano Zampini } 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE)); 145ee1cef2cSJed Brown PetscFunctionReturn(0); 146ee1cef2cSJed Brown } 147ee1cef2cSJed Brown 148ee1cef2cSJed Brown static PetscErrorCode MatDestroy_SubMatrix(Mat N) 149ee1cef2cSJed Brown { 15054e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 151ee1cef2cSJed Brown 152ee1cef2cSJed Brown PetscFunctionBegin; 1535f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&Na->isrow)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&Na->iscol)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->lwork)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->rwork)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->lwork2)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->rwork2)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&Na->lrestrict)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&Na->rprolong)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Na->A)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(N->data)); 163ee1cef2cSJed Brown PetscFunctionReturn(0); 164ee1cef2cSJed Brown } 165ee1cef2cSJed Brown 166ee1cef2cSJed Brown /*@ 16754e05e6cSHong Zhang MatCreateSubMatrixVirtual - Creates a virtual matrix that acts as a submatrix 168ee1cef2cSJed Brown 169ee1cef2cSJed Brown Collective on Mat 170ee1cef2cSJed Brown 171ee1cef2cSJed Brown Input Parameters: 172ee1cef2cSJed Brown + A - matrix that we will extract a submatrix of 173ee1cef2cSJed Brown . isrow - rows to be present in the submatrix 174ee1cef2cSJed Brown - iscol - columns to be present in the submatrix 175ee1cef2cSJed Brown 176ee1cef2cSJed Brown Output Parameters: 177ee1cef2cSJed Brown . newmat - new matrix 178ee1cef2cSJed Brown 179ee1cef2cSJed Brown Level: developer 180ee1cef2cSJed Brown 181ee1cef2cSJed Brown Notes: 1827dae84e0SHong Zhang Most will use MatCreateSubMatrix which provides a more efficient representation if it is available. 183ee1cef2cSJed Brown 18454e05e6cSHong Zhang .seealso: MatCreateSubMatrix(), MatSubMatrixVirtualUpdate() 185ee1cef2cSJed Brown @*/ 18654e05e6cSHong Zhang PetscErrorCode MatCreateSubMatrixVirtual(Mat A,IS isrow,IS iscol,Mat *newmat) 187ee1cef2cSJed Brown { 188ee1cef2cSJed Brown Vec left,right; 189ee1cef2cSJed Brown PetscInt m,n; 190ee1cef2cSJed Brown Mat N; 19154e05e6cSHong Zhang Mat_SubVirtual *Na; 192ee1cef2cSJed Brown 193ee1cef2cSJed Brown PetscFunctionBegin; 1940700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1950700a824SBarry Smith PetscValidHeaderSpecific(isrow,IS_CLASSID,2); 1960700a824SBarry Smith PetscValidHeaderSpecific(iscol,IS_CLASSID,3); 197ee1cef2cSJed Brown PetscValidPointer(newmat,4); 198f4259b30SLisandro Dalcin *newmat = NULL; 199ee1cef2cSJed Brown 2005f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),&N)); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrow,&m)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscol,&n)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(N,m,n,PETSC_DETERMINE,PETSC_DETERMINE)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)N,MATSUBMATRIX)); 205ee1cef2cSJed Brown 2065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(N,&Na)); 207ee1cef2cSJed Brown N->data = (void*)Na; 20811c5f74dSStefano Zampini 2095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)isrow)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)iscol)); 211ee1cef2cSJed Brown Na->isrow = isrow; 212ee1cef2cSJed Brown Na->iscol = iscol; 21311c5f74dSStefano Zampini 2145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(N->defaultvectype)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(A->defaultvectype,&N->defaultvectype)); 21611c5f74dSStefano Zampini /* Do not use MatConvert directly since MatShell has a duplicate operation which does not increase 21711c5f74dSStefano Zampini the reference count of the context. This is a problem if A is already of type MATSHELL */ 2185f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvertFrom_Shell(A,MATSHELL,MAT_INITIAL_MATRIX,&Na->A)); 219ee1cef2cSJed Brown 220ee1cef2cSJed Brown N->ops->destroy = MatDestroy_SubMatrix; 221ee1cef2cSJed Brown N->ops->mult = MatMult_SubMatrix; 222ee1cef2cSJed Brown N->ops->multadd = MatMultAdd_SubMatrix; 223ee1cef2cSJed Brown N->ops->multtranspose = MatMultTranspose_SubMatrix; 224ee1cef2cSJed Brown N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix; 225ee1cef2cSJed Brown N->ops->scale = MatScale_SubMatrix; 226ee1cef2cSJed Brown N->ops->diagonalscale = MatDiagonalScale_SubMatrix; 227bddd1ffdSAlp Dener N->ops->shift = MatShift_SubMatrix; 22811c5f74dSStefano Zampini N->ops->convert = MatConvert_Shell; 22911c5f74dSStefano Zampini N->ops->getdiagonal = MatGetDiagonal_SubMatrix; 230ee1cef2cSJed Brown 2315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSizesFromMats(N,A,A)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetUp(N->rmap)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutSetUp(N->cmap)); 234ee1cef2cSJed Brown 2355f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&Na->rwork,&Na->lwork)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(N,&right,&left)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(Na->lwork,isrow,left,NULL,&Na->lrestrict)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(right,NULL,Na->rwork,iscol,&Na->rprolong)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&left)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&right)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(N)); 24226fbe8dcSKarl Rupp 24311c5f74dSStefano Zampini N->assembled = PETSC_TRUE; 244ee1cef2cSJed Brown *newmat = N; 245ee1cef2cSJed Brown PetscFunctionReturn(0); 246ee1cef2cSJed Brown } 247ee1cef2cSJed Brown 248ee1cef2cSJed Brown /*@ 24954e05e6cSHong Zhang MatSubMatrixVirtualUpdate - Updates a submatrix 250ee1cef2cSJed Brown 251ee1cef2cSJed Brown Collective on Mat 252ee1cef2cSJed Brown 253ee1cef2cSJed Brown Input Parameters: 254ee1cef2cSJed Brown + N - submatrix to update 255ee1cef2cSJed Brown . A - full matrix in the submatrix 256ee1cef2cSJed Brown . isrow - rows in the update (same as the first time the submatrix was created) 257ee1cef2cSJed Brown - iscol - columns in the update (same as the first time the submatrix was created) 258ee1cef2cSJed Brown 259ee1cef2cSJed Brown Level: developer 260ee1cef2cSJed Brown 261ee1cef2cSJed Brown Notes: 2627dae84e0SHong Zhang Most will use MatCreateSubMatrix which provides a more efficient representation if it is available. 263ee1cef2cSJed Brown 26454e05e6cSHong Zhang .seealso: MatCreateSubMatrixVirtual() 265ee1cef2cSJed Brown @*/ 26654e05e6cSHong Zhang PetscErrorCode MatSubMatrixVirtualUpdate(Mat N,Mat A,IS isrow,IS iscol) 267ee1cef2cSJed Brown { 268ace3abfcSBarry Smith PetscBool flg; 26954e05e6cSHong Zhang Mat_SubVirtual *Na; 270ee1cef2cSJed Brown 271ee1cef2cSJed Brown PetscFunctionBegin; 2720700a824SBarry Smith PetscValidHeaderSpecific(N,MAT_CLASSID,1); 2730700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,2); 2740700a824SBarry Smith PetscValidHeaderSpecific(isrow,IS_CLASSID,3); 2750700a824SBarry Smith PetscValidHeaderSpecific(iscol,IS_CLASSID,4); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)N,MATSUBMATRIX,&flg)); 277*28b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix has wrong type"); 278ee1cef2cSJed Brown 27954e05e6cSHong Zhang Na = (Mat_SubVirtual*)N->data; 2805f80ce2aSJacob Faibussowitsch CHKERRQ(ISEqual(isrow,Na->isrow,&flg)); 281*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different row indices"); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(ISEqual(iscol,Na->iscol,&flg)); 283*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different column indices"); 284ee1cef2cSJed Brown 2855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(N->defaultvectype)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(A->defaultvectype,&N->defaultvectype)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Na->A)); 28811c5f74dSStefano Zampini /* Do not use MatConvert directly since MatShell has a duplicate operation which does not increase 28911c5f74dSStefano Zampini the reference count of the context. This is a problem if A is already of type MATSHELL */ 2905f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvertFrom_Shell(A,MATSHELL,MAT_INITIAL_MATRIX,&Na->A)); 291ee1cef2cSJed Brown PetscFunctionReturn(0); 292ee1cef2cSJed Brown } 293