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 129371c9d4SSatish Balay static PetscErrorCode MatScale_SubMatrix(Mat N, PetscScalar a) { 13bddd1ffdSAlp Dener Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 14bddd1ffdSAlp Dener 15bddd1ffdSAlp Dener PetscFunctionBegin; 169566063dSJacob Faibussowitsch PetscCall(MatScale(Na->A, a)); 17bddd1ffdSAlp Dener PetscFunctionReturn(0); 18bddd1ffdSAlp Dener } 19bddd1ffdSAlp Dener 209371c9d4SSatish Balay static PetscErrorCode MatShift_SubMatrix(Mat N, PetscScalar a) { 21bddd1ffdSAlp Dener Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 22bddd1ffdSAlp Dener 23bddd1ffdSAlp Dener PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(MatShift(Na->A, a)); 25ee1cef2cSJed Brown PetscFunctionReturn(0); 26ee1cef2cSJed Brown } 27ee1cef2cSJed Brown 289371c9d4SSatish Balay static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N, Vec left, Vec right) { 2954e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 30ee1cef2cSJed Brown 31ee1cef2cSJed Brown PetscFunctionBegin; 32ee1cef2cSJed Brown if (right) { 339566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork)); 349566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, right, Na->rwork, INSERT_VALUES, SCATTER_FORWARD)); 359566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, right, Na->rwork, INSERT_VALUES, SCATTER_FORWARD)); 36ee1cef2cSJed Brown } 3711c5f74dSStefano Zampini if (left) { 389566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork)); 399566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, left, Na->lwork, INSERT_VALUES, SCATTER_REVERSE)); 409566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, left, Na->lwork, INSERT_VALUES, SCATTER_REVERSE)); 41ee1cef2cSJed Brown } 429566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(Na->A, left ? Na->lwork : NULL, right ? Na->rwork : NULL)); 4311c5f74dSStefano Zampini PetscFunctionReturn(0); 4411c5f74dSStefano Zampini } 4511c5f74dSStefano Zampini 469371c9d4SSatish Balay static PetscErrorCode MatGetDiagonal_SubMatrix(Mat N, Vec d) { 4711c5f74dSStefano Zampini Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 4811c5f74dSStefano Zampini 4911c5f74dSStefano Zampini PetscFunctionBegin; 509566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(Na->A, Na->rwork)); 519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, Na->rwork, d, INSERT_VALUES, SCATTER_REVERSE)); 529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, Na->rwork, d, INSERT_VALUES, SCATTER_REVERSE)); 53ee1cef2cSJed Brown PetscFunctionReturn(0); 54ee1cef2cSJed Brown } 55ee1cef2cSJed Brown 569371c9d4SSatish Balay static PetscErrorCode MatMult_SubMatrix(Mat N, Vec x, Vec y) { 5754e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 58ee1cef2cSJed Brown 59ee1cef2cSJed Brown PetscFunctionBegin; 609566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork)); 619566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, x, Na->rwork, INSERT_VALUES, SCATTER_FORWARD)); 629566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, x, Na->rwork, INSERT_VALUES, SCATTER_FORWARD)); 639566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, Na->rwork, Na->lwork)); 649566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, Na->lwork, y, INSERT_VALUES, SCATTER_FORWARD)); 659566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, Na->lwork, y, INSERT_VALUES, SCATTER_FORWARD)); 66ee1cef2cSJed Brown PetscFunctionReturn(0); 67ee1cef2cSJed Brown } 68ee1cef2cSJed Brown 699371c9d4SSatish Balay static PetscErrorCode MatMultAdd_SubMatrix(Mat N, Vec v1, Vec v2, Vec v3) { 7054e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 71ee1cef2cSJed Brown 72ee1cef2cSJed Brown PetscFunctionBegin; 739566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork)); 749566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, v1, Na->rwork, INSERT_VALUES, SCATTER_FORWARD)); 759566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, v1, Na->rwork, INSERT_VALUES, SCATTER_FORWARD)); 7611c5f74dSStefano Zampini if (v1 == v2) { 779566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Na->A, Na->rwork, Na->rwork, Na->lwork)); 7811c5f74dSStefano Zampini } else if (v2 == v3) { 799566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork)); 809566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, v2, Na->lwork, INSERT_VALUES, SCATTER_REVERSE)); 819566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, v2, Na->lwork, INSERT_VALUES, SCATTER_REVERSE)); 829566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Na->A, Na->rwork, Na->lwork, Na->lwork)); 831e7f0bbdSJed Brown } else { 8411c5f74dSStefano Zampini if (!Na->lwork2) { 859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->lwork, &Na->lwork2)); 8611c5f74dSStefano Zampini } else { 879566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork2)); 8811c5f74dSStefano Zampini } 899566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, v2, Na->lwork2, INSERT_VALUES, SCATTER_REVERSE)); 909566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, v2, Na->lwork2, INSERT_VALUES, SCATTER_REVERSE)); 919566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Na->A, Na->rwork, Na->lwork2, Na->lwork)); 9211c5f74dSStefano Zampini } 939566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, Na->lwork, v3, INSERT_VALUES, SCATTER_FORWARD)); 949566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, Na->lwork, v3, INSERT_VALUES, SCATTER_FORWARD)); 95ee1cef2cSJed Brown PetscFunctionReturn(0); 96ee1cef2cSJed Brown } 97ee1cef2cSJed Brown 989371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_SubMatrix(Mat N, Vec x, Vec y) { 9954e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 100ee1cef2cSJed Brown 101ee1cef2cSJed Brown PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork)); 1039566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, x, Na->lwork, INSERT_VALUES, SCATTER_REVERSE)); 1049566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, x, Na->lwork, INSERT_VALUES, SCATTER_REVERSE)); 1059566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A, Na->lwork, Na->rwork)); 1069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, Na->rwork, y, INSERT_VALUES, SCATTER_REVERSE)); 1079566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, Na->rwork, y, INSERT_VALUES, SCATTER_REVERSE)); 108ee1cef2cSJed Brown PetscFunctionReturn(0); 109ee1cef2cSJed Brown } 110ee1cef2cSJed Brown 1119371c9d4SSatish Balay static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N, Vec v1, Vec v2, Vec v3) { 11254e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 113ee1cef2cSJed Brown 114ee1cef2cSJed Brown PetscFunctionBegin; 1159566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork)); 1169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, v1, Na->lwork, INSERT_VALUES, SCATTER_REVERSE)); 1179566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, v1, Na->lwork, INSERT_VALUES, SCATTER_REVERSE)); 11811c5f74dSStefano Zampini if (v1 == v2) { 1199566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(Na->A, Na->lwork, Na->lwork, Na->rwork)); 12011c5f74dSStefano Zampini } else if (v2 == v3) { 1219566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork)); 1229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, v2, Na->rwork, INSERT_VALUES, SCATTER_FORWARD)); 1239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, v2, Na->rwork, INSERT_VALUES, SCATTER_FORWARD)); 1249566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(Na->A, Na->lwork, Na->rwork, Na->rwork)); 1251e7f0bbdSJed Brown } else { 12611c5f74dSStefano Zampini if (!Na->rwork2) { 1279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->rwork, &Na->rwork2)); 12811c5f74dSStefano Zampini } else { 1299566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork2)); 13011c5f74dSStefano Zampini } 1319566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, v2, Na->rwork2, INSERT_VALUES, SCATTER_FORWARD)); 1329566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, v2, Na->rwork2, INSERT_VALUES, SCATTER_FORWARD)); 1339566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(Na->A, Na->lwork, Na->rwork2, Na->rwork)); 13411c5f74dSStefano Zampini } 1359566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, Na->rwork, v3, INSERT_VALUES, SCATTER_REVERSE)); 1369566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, Na->rwork, v3, INSERT_VALUES, SCATTER_REVERSE)); 137ee1cef2cSJed Brown PetscFunctionReturn(0); 138ee1cef2cSJed Brown } 139ee1cef2cSJed Brown 1409371c9d4SSatish Balay static PetscErrorCode MatDestroy_SubMatrix(Mat N) { 14154e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data; 142ee1cef2cSJed Brown 143ee1cef2cSJed Brown PetscFunctionBegin; 1449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&Na->isrow)); 1459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&Na->iscol)); 1469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->lwork)); 1479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rwork)); 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->lwork2)); 1499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rwork2)); 1509566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&Na->lrestrict)); 1519566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&Na->rprolong)); 1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 1539566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 154ee1cef2cSJed Brown PetscFunctionReturn(0); 155ee1cef2cSJed Brown } 156ee1cef2cSJed Brown 157ee1cef2cSJed Brown /*@ 15811a5261eSBarry Smith MatCreateSubMatrixVirtual - Creates a virtual matrix `MATSUBMATRIX` that acts as a submatrix 159ee1cef2cSJed Brown 16011a5261eSBarry Smith Collective on A 161ee1cef2cSJed Brown 162ee1cef2cSJed Brown Input Parameters: 163ee1cef2cSJed Brown + A - matrix that we will extract a submatrix of 164ee1cef2cSJed Brown . isrow - rows to be present in the submatrix 165ee1cef2cSJed Brown - iscol - columns to be present in the submatrix 166ee1cef2cSJed Brown 167ee1cef2cSJed Brown Output Parameters: 168ee1cef2cSJed Brown . newmat - new matrix 169ee1cef2cSJed Brown 170ee1cef2cSJed Brown Level: developer 171ee1cef2cSJed Brown 17211a5261eSBarry Smith Note: 17311a5261eSBarry Smith Most will use `MatCreateSubMatrix()` which provides a more efficient representation if it is available. 174ee1cef2cSJed Brown 17511a5261eSBarry Smith Developer Note: 17611a5261eSBarry Smith The `MatType` is `MATSUBMATRIX` but the routines associated have `SubMatrixVirtual` in them, the `MatType` should likely be changed 17711a5261eSBarry Smith 17811a5261eSBarry Smith .seealso: `MATSUBMATRIX`, `MATLOCALREF`, `MatCreateLocalRef()`, `MatCreateSubMatrix()`, `MatSubMatrixVirtualUpdate()` 179ee1cef2cSJed Brown @*/ 1809371c9d4SSatish Balay PetscErrorCode MatCreateSubMatrixVirtual(Mat A, IS isrow, IS iscol, Mat *newmat) { 181ee1cef2cSJed Brown Vec left, right; 182ee1cef2cSJed Brown PetscInt m, n; 183ee1cef2cSJed Brown Mat N; 18454e05e6cSHong Zhang Mat_SubVirtual *Na; 185ee1cef2cSJed Brown 186ee1cef2cSJed Brown PetscFunctionBegin; 1870700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 1880700a824SBarry Smith PetscValidHeaderSpecific(isrow, IS_CLASSID, 2); 1890700a824SBarry Smith PetscValidHeaderSpecific(iscol, IS_CLASSID, 3); 190ee1cef2cSJed Brown PetscValidPointer(newmat, 4); 191f4259b30SLisandro Dalcin *newmat = NULL; 192ee1cef2cSJed Brown 1939566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &N)); 1949566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow, &m)); 1959566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscol, &n)); 1969566063dSJacob Faibussowitsch PetscCall(MatSetSizes(N, m, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1979566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)N, MATSUBMATRIX)); 198ee1cef2cSJed Brown 199*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&Na)); 200ee1cef2cSJed Brown N->data = (void *)Na; 20111c5f74dSStefano Zampini 2029566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 2039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 204ee1cef2cSJed Brown Na->isrow = isrow; 205ee1cef2cSJed Brown Na->iscol = iscol; 20611c5f74dSStefano Zampini 2079566063dSJacob Faibussowitsch PetscCall(PetscFree(N->defaultvectype)); 2089566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(A->defaultvectype, &N->defaultvectype)); 20911c5f74dSStefano Zampini /* Do not use MatConvert directly since MatShell has a duplicate operation which does not increase 21011c5f74dSStefano Zampini the reference count of the context. This is a problem if A is already of type MATSHELL */ 2119566063dSJacob Faibussowitsch PetscCall(MatConvertFrom_Shell(A, MATSHELL, MAT_INITIAL_MATRIX, &Na->A)); 212ee1cef2cSJed Brown 213ee1cef2cSJed Brown N->ops->destroy = MatDestroy_SubMatrix; 214ee1cef2cSJed Brown N->ops->mult = MatMult_SubMatrix; 215ee1cef2cSJed Brown N->ops->multadd = MatMultAdd_SubMatrix; 216ee1cef2cSJed Brown N->ops->multtranspose = MatMultTranspose_SubMatrix; 217ee1cef2cSJed Brown N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix; 218ee1cef2cSJed Brown N->ops->scale = MatScale_SubMatrix; 219ee1cef2cSJed Brown N->ops->diagonalscale = MatDiagonalScale_SubMatrix; 220bddd1ffdSAlp Dener N->ops->shift = MatShift_SubMatrix; 22111c5f74dSStefano Zampini N->ops->convert = MatConvert_Shell; 22211c5f74dSStefano Zampini N->ops->getdiagonal = MatGetDiagonal_SubMatrix; 223ee1cef2cSJed Brown 2249566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(N, A, A)); 2259566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(N->rmap)); 2269566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(N->cmap)); 227ee1cef2cSJed Brown 2289566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Na->rwork, &Na->lwork)); 2299566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(N, &right, &left)); 2309566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(Na->lwork, isrow, left, NULL, &Na->lrestrict)); 2319566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(right, NULL, Na->rwork, iscol, &Na->rprolong)); 2329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left)); 2339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right)); 2349566063dSJacob Faibussowitsch PetscCall(MatSetUp(N)); 23526fbe8dcSKarl Rupp 23611c5f74dSStefano Zampini N->assembled = PETSC_TRUE; 237ee1cef2cSJed Brown *newmat = N; 238ee1cef2cSJed Brown PetscFunctionReturn(0); 239ee1cef2cSJed Brown } 240ee1cef2cSJed Brown 24111a5261eSBarry Smith /*MC 24211a5261eSBarry Smith MATSUBMATRIX - "submatrix" - A matrix type that represents a virtual submatrix of a matrix 243ee1cef2cSJed Brown 24411a5261eSBarry Smith Level: advanced 24511a5261eSBarry Smith 24611a5261eSBarry Smith Developer Note: 24711a5261eSBarry Smith This should be named `MATSUBMATRIXVIRTUAL` 24811a5261eSBarry Smith 24911a5261eSBarry Smith .seealso: `Mat`, `MatCreateSubMatrixVirtual()`, `MatCreateSubMatrixVirtual()`, `MatCreateSubMatrix()` 25011a5261eSBarry Smith M*/ 25111a5261eSBarry Smith 25211a5261eSBarry Smith /*@ 25311a5261eSBarry Smith MatSubMatrixVirtualUpdate - Updates a `MATSUBMATRIX` virtual submatrix 25411a5261eSBarry Smith 25511a5261eSBarry Smith Collective on N 256ee1cef2cSJed Brown 257ee1cef2cSJed Brown Input Parameters: 258ee1cef2cSJed Brown + N - submatrix to update 259ee1cef2cSJed Brown . A - full matrix in the submatrix 260ee1cef2cSJed Brown . isrow - rows in the update (same as the first time the submatrix was created) 261ee1cef2cSJed Brown - iscol - columns in the update (same as the first time the submatrix was created) 262ee1cef2cSJed Brown 263ee1cef2cSJed Brown Level: developer 264ee1cef2cSJed Brown 26511a5261eSBarry Smith Note: 26611a5261eSBarry Smith Most will use `MatCreateSubMatrix()` which provides a more efficient representation if it is available. 267ee1cef2cSJed Brown 26811a5261eSBarry Smith .seealso: MATSUBMATRIX`, `MatCreateSubMatrixVirtual()` 269ee1cef2cSJed Brown @*/ 2709371c9d4SSatish Balay PetscErrorCode MatSubMatrixVirtualUpdate(Mat N, Mat A, IS isrow, IS iscol) { 271ace3abfcSBarry Smith PetscBool flg; 27254e05e6cSHong Zhang Mat_SubVirtual *Na; 273ee1cef2cSJed Brown 274ee1cef2cSJed Brown PetscFunctionBegin; 2750700a824SBarry Smith PetscValidHeaderSpecific(N, MAT_CLASSID, 1); 2760700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 2); 2770700a824SBarry Smith PetscValidHeaderSpecific(isrow, IS_CLASSID, 3); 2780700a824SBarry Smith PetscValidHeaderSpecific(iscol, IS_CLASSID, 4); 2799566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)N, MATSUBMATRIX, &flg)); 28028b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix has wrong type"); 281ee1cef2cSJed Brown 28254e05e6cSHong Zhang Na = (Mat_SubVirtual *)N->data; 2839566063dSJacob Faibussowitsch PetscCall(ISEqual(isrow, Na->isrow, &flg)); 28428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot update submatrix with different row indices"); 2859566063dSJacob Faibussowitsch PetscCall(ISEqual(iscol, Na->iscol, &flg)); 28628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot update submatrix with different column indices"); 287ee1cef2cSJed Brown 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(N->defaultvectype)); 2899566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(A->defaultvectype, &N->defaultvectype)); 2909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 29111c5f74dSStefano Zampini /* Do not use MatConvert directly since MatShell has a duplicate operation which does not increase 29211c5f74dSStefano Zampini the reference count of the context. This is a problem if A is already of type MATSHELL */ 2939566063dSJacob Faibussowitsch PetscCall(MatConvertFrom_Shell(A, MATSHELL, MAT_INITIAL_MATRIX, &Na->A)); 294ee1cef2cSJed Brown PetscFunctionReturn(0); 295ee1cef2cSJed Brown } 296