xref: /petsc/src/mat/impls/submat/submat.c (revision db7814771ca77b190574494e87b584e981451db0)
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;
179566063dSJacob Faibussowitsch   PetscCall(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;
269566063dSJacob Faibussowitsch   PetscCall(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) {
369566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Na->rwork));
379566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(Na->rprolong,right,Na->rwork,INSERT_VALUES,SCATTER_FORWARD));
389566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(Na->rprolong,right,Na->rwork,INSERT_VALUES,SCATTER_FORWARD));
39ee1cef2cSJed Brown   }
4011c5f74dSStefano Zampini   if (left) {
419566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Na->lwork));
429566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(Na->lrestrict,left,Na->lwork,INSERT_VALUES,SCATTER_REVERSE));
439566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(Na->lrestrict,left,Na->lwork,INSERT_VALUES,SCATTER_REVERSE));
44ee1cef2cSJed Brown   }
459566063dSJacob Faibussowitsch   PetscCall(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;
549566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(Na->A,Na->rwork));
559566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->rprolong,Na->rwork,d,INSERT_VALUES,SCATTER_REVERSE));
569566063dSJacob Faibussowitsch   PetscCall(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;
659566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Na->rwork));
669566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->rprolong,x,Na->rwork,INSERT_VALUES,SCATTER_FORWARD));
679566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(Na->rprolong,x,Na->rwork,INSERT_VALUES,SCATTER_FORWARD));
689566063dSJacob Faibussowitsch   PetscCall(MatMult(Na->A,Na->rwork,Na->lwork));
699566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD));
709566063dSJacob Faibussowitsch   PetscCall(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;
799566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Na->rwork));
809566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->rprolong,v1,Na->rwork,INSERT_VALUES,SCATTER_FORWARD));
819566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(Na->rprolong,v1,Na->rwork,INSERT_VALUES,SCATTER_FORWARD));
8211c5f74dSStefano Zampini   if (v1 == v2) {
839566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(Na->A,Na->rwork,Na->rwork,Na->lwork));
8411c5f74dSStefano Zampini   } else if (v2 == v3) {
859566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Na->lwork));
869566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(Na->lrestrict,v2,Na->lwork,INSERT_VALUES,SCATTER_REVERSE));
879566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(Na->lrestrict,v2,Na->lwork,INSERT_VALUES,SCATTER_REVERSE));
889566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(Na->A,Na->rwork,Na->lwork,Na->lwork));
891e7f0bbdSJed Brown   } else {
9011c5f74dSStefano Zampini     if (!Na->lwork2) {
919566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(Na->lwork,&Na->lwork2));
9211c5f74dSStefano Zampini     } else {
939566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Na->lwork2));
9411c5f74dSStefano Zampini     }
959566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(Na->lrestrict,v2,Na->lwork2,INSERT_VALUES,SCATTER_REVERSE));
969566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(Na->lrestrict,v2,Na->lwork2,INSERT_VALUES,SCATTER_REVERSE));
979566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(Na->A,Na->rwork,Na->lwork2,Na->lwork));
9811c5f74dSStefano Zampini   }
999566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD));
1009566063dSJacob Faibussowitsch   PetscCall(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;
1099566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Na->lwork));
1109566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->lrestrict,x,Na->lwork,INSERT_VALUES,SCATTER_REVERSE));
1119566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(Na->lrestrict,x,Na->lwork,INSERT_VALUES,SCATTER_REVERSE));
1129566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(Na->A,Na->lwork,Na->rwork));
1139566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE));
1149566063dSJacob Faibussowitsch   PetscCall(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;
1239566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Na->lwork));
1249566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->lrestrict,v1,Na->lwork,INSERT_VALUES,SCATTER_REVERSE));
1259566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(Na->lrestrict,v1,Na->lwork,INSERT_VALUES,SCATTER_REVERSE));
12611c5f74dSStefano Zampini   if (v1 == v2) {
1279566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(Na->A,Na->lwork,Na->lwork,Na->rwork));
12811c5f74dSStefano Zampini   } else if (v2 == v3) {
1299566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(Na->rwork));
1309566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(Na->rprolong,v2,Na->rwork,INSERT_VALUES,SCATTER_FORWARD));
1319566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(Na->rprolong,v2,Na->rwork,INSERT_VALUES,SCATTER_FORWARD));
1329566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(Na->A,Na->lwork,Na->rwork,Na->rwork));
1331e7f0bbdSJed Brown   } else {
13411c5f74dSStefano Zampini     if (!Na->rwork2) {
1359566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(Na->rwork,&Na->rwork2));
13611c5f74dSStefano Zampini     } else {
1379566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(Na->rwork2));
13811c5f74dSStefano Zampini     }
1399566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(Na->rprolong,v2,Na->rwork2,INSERT_VALUES,SCATTER_FORWARD));
1409566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(Na->rprolong,v2,Na->rwork2,INSERT_VALUES,SCATTER_FORWARD));
1419566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(Na->A,Na->lwork,Na->rwork2,Na->rwork));
14211c5f74dSStefano Zampini   }
1439566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE));
1449566063dSJacob Faibussowitsch   PetscCall(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;
1539566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Na->isrow));
1549566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&Na->iscol));
1559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->lwork));
1569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->rwork));
1579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->lwork2));
1589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Na->rwork2));
1599566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&Na->lrestrict));
1609566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&Na->rprolong));
1619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Na->A));
1629566063dSJacob Faibussowitsch   PetscCall(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 
184*db781477SPatrick Sanan .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 
2009566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&N));
2019566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrow,&m));
2029566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscol,&n));
2039566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(N,m,n,PETSC_DETERMINE,PETSC_DETERMINE));
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)N,MATSUBMATRIX));
205ee1cef2cSJed Brown 
2069566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(N,&Na));
207ee1cef2cSJed Brown   N->data   = (void*)Na;
20811c5f74dSStefano Zampini 
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)isrow));
2109566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)iscol));
211ee1cef2cSJed Brown   Na->isrow = isrow;
212ee1cef2cSJed Brown   Na->iscol = iscol;
21311c5f74dSStefano Zampini 
2149566063dSJacob Faibussowitsch   PetscCall(PetscFree(N->defaultvectype));
2159566063dSJacob Faibussowitsch   PetscCall(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 */
2189566063dSJacob Faibussowitsch   PetscCall(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 
2319566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(N,A,A));
2329566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(N->rmap));
2339566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(N->cmap));
234ee1cef2cSJed Brown 
2359566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&Na->rwork,&Na->lwork));
2369566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(N,&right,&left));
2379566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(Na->lwork,isrow,left,NULL,&Na->lrestrict));
2389566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(right,NULL,Na->rwork,iscol,&Na->rprolong));
2399566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&left));
2409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&right));
2419566063dSJacob Faibussowitsch   PetscCall(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 
264*db781477SPatrick Sanan .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);
2769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)N,MATSUBMATRIX,&flg));
27728b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix has wrong type");
278ee1cef2cSJed Brown 
27954e05e6cSHong Zhang   Na   = (Mat_SubVirtual*)N->data;
2809566063dSJacob Faibussowitsch   PetscCall(ISEqual(isrow,Na->isrow,&flg));
28128b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different row indices");
2829566063dSJacob Faibussowitsch   PetscCall(ISEqual(iscol,Na->iscol,&flg));
28328b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different column indices");
284ee1cef2cSJed Brown 
2859566063dSJacob Faibussowitsch   PetscCall(PetscFree(N->defaultvectype));
2869566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(A->defaultvectype,&N->defaultvectype));
2879566063dSJacob Faibussowitsch   PetscCall(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 */
2909566063dSJacob Faibussowitsch   PetscCall(MatConvertFrom_Shell(A,MATSHELL,MAT_INITIAL_MATRIX,&Na->A));
291ee1cef2cSJed Brown   PetscFunctionReturn(0);
292ee1cef2cSJed Brown }
293