1af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
2ee1cef2cSJed Brown
3ee1cef2cSJed Brown typedef struct {
4ee1cef2cSJed Brown IS isrow, iscol; /* rows and columns in submatrix, only used to check consistency */
5ee1cef2cSJed Brown Vec lwork, rwork; /* work vectors inside the scatters */
611c5f74dSStefano Zampini Vec lwork2, rwork2; /* work vectors inside the scatters */
7ee1cef2cSJed Brown VecScatter lrestrict, rprolong;
8ee1cef2cSJed Brown Mat A;
954e05e6cSHong Zhang } Mat_SubVirtual;
10ee1cef2cSJed Brown
MatScale_SubMatrix(Mat N,PetscScalar a)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_SubMatrix(Mat N, PetscScalar a)
12d71ae5a4SJacob Faibussowitsch {
13bddd1ffdSAlp Dener Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
14bddd1ffdSAlp Dener
15bddd1ffdSAlp Dener PetscFunctionBegin;
169566063dSJacob Faibussowitsch PetscCall(MatScale(Na->A, a));
173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
18bddd1ffdSAlp Dener }
19bddd1ffdSAlp Dener
MatShift_SubMatrix(Mat N,PetscScalar a)20d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatShift_SubMatrix(Mat N, PetscScalar a)
21d71ae5a4SJacob Faibussowitsch {
22bddd1ffdSAlp Dener Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
23bddd1ffdSAlp Dener
24bddd1ffdSAlp Dener PetscFunctionBegin;
259566063dSJacob Faibussowitsch PetscCall(MatShift(Na->A, a));
263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27ee1cef2cSJed Brown }
28ee1cef2cSJed Brown
MatDiagonalScale_SubMatrix(Mat N,Vec left,Vec right)29d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N, Vec left, Vec right)
30d71ae5a4SJacob Faibussowitsch {
3154e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
32ee1cef2cSJed Brown
33ee1cef2cSJed Brown PetscFunctionBegin;
34ee1cef2cSJed Brown if (right) {
359566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork));
369566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, right, Na->rwork, INSERT_VALUES, SCATTER_FORWARD));
379566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, right, Na->rwork, INSERT_VALUES, SCATTER_FORWARD));
38ee1cef2cSJed Brown }
3911c5f74dSStefano Zampini if (left) {
409566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork));
419566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, left, Na->lwork, INSERT_VALUES, SCATTER_REVERSE));
429566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, left, Na->lwork, INSERT_VALUES, SCATTER_REVERSE));
43ee1cef2cSJed Brown }
449566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(Na->A, left ? Na->lwork : NULL, right ? Na->rwork : NULL));
453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4611c5f74dSStefano Zampini }
4711c5f74dSStefano Zampini
MatGetDiagonal_SubMatrix(Mat N,Vec d)48d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_SubMatrix(Mat N, Vec d)
49d71ae5a4SJacob Faibussowitsch {
5011c5f74dSStefano Zampini Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
5111c5f74dSStefano Zampini
5211c5f74dSStefano Zampini PetscFunctionBegin;
539566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(Na->A, Na->rwork));
549566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, Na->rwork, d, INSERT_VALUES, SCATTER_REVERSE));
559566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, Na->rwork, d, INSERT_VALUES, SCATTER_REVERSE));
563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
57ee1cef2cSJed Brown }
58ee1cef2cSJed Brown
MatMult_SubMatrix(Mat N,Vec x,Vec y)59d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_SubMatrix(Mat N, Vec x, Vec y)
60d71ae5a4SJacob Faibussowitsch {
6154e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
62ee1cef2cSJed Brown
63ee1cef2cSJed Brown PetscFunctionBegin;
649566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork));
659566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, x, Na->rwork, INSERT_VALUES, SCATTER_FORWARD));
669566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, x, Na->rwork, INSERT_VALUES, SCATTER_FORWARD));
679566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A, Na->rwork, Na->lwork));
689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, Na->lwork, y, INSERT_VALUES, SCATTER_FORWARD));
699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, Na->lwork, y, INSERT_VALUES, SCATTER_FORWARD));
703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
71ee1cef2cSJed Brown }
72ee1cef2cSJed Brown
MatMultAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3)73d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_SubMatrix(Mat N, Vec v1, Vec v2, Vec v3)
74d71ae5a4SJacob Faibussowitsch {
7554e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
76ee1cef2cSJed Brown
77ee1cef2cSJed Brown PetscFunctionBegin;
789566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork));
799566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, v1, Na->rwork, INSERT_VALUES, SCATTER_FORWARD));
809566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, v1, Na->rwork, INSERT_VALUES, SCATTER_FORWARD));
8111c5f74dSStefano Zampini if (v1 == v2) {
829566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Na->A, Na->rwork, Na->rwork, Na->lwork));
8311c5f74dSStefano Zampini } else if (v2 == v3) {
849566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork));
859566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, v2, Na->lwork, INSERT_VALUES, SCATTER_REVERSE));
869566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, v2, Na->lwork, INSERT_VALUES, SCATTER_REVERSE));
879566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Na->A, Na->rwork, Na->lwork, Na->lwork));
881e7f0bbdSJed Brown } else {
8911c5f74dSStefano Zampini if (!Na->lwork2) {
909566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->lwork, &Na->lwork2));
9111c5f74dSStefano Zampini } else {
929566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork2));
9311c5f74dSStefano Zampini }
949566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, v2, Na->lwork2, INSERT_VALUES, SCATTER_REVERSE));
959566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, v2, Na->lwork2, INSERT_VALUES, SCATTER_REVERSE));
969566063dSJacob Faibussowitsch PetscCall(MatMultAdd(Na->A, Na->rwork, Na->lwork2, Na->lwork));
9711c5f74dSStefano Zampini }
989566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, Na->lwork, v3, INSERT_VALUES, SCATTER_FORWARD));
999566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, Na->lwork, v3, INSERT_VALUES, SCATTER_FORWARD));
1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
101ee1cef2cSJed Brown }
102ee1cef2cSJed Brown
MatMultTranspose_SubMatrix(Mat N,Vec x,Vec y)103d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_SubMatrix(Mat N, Vec x, Vec y)
104d71ae5a4SJacob Faibussowitsch {
10554e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
106ee1cef2cSJed Brown
107ee1cef2cSJed Brown PetscFunctionBegin;
1089566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork));
1099566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, x, Na->lwork, INSERT_VALUES, SCATTER_REVERSE));
1109566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, x, Na->lwork, INSERT_VALUES, SCATTER_REVERSE));
1119566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A, Na->lwork, Na->rwork));
1129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, Na->rwork, y, INSERT_VALUES, SCATTER_REVERSE));
1139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, Na->rwork, y, INSERT_VALUES, SCATTER_REVERSE));
1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
115ee1cef2cSJed Brown }
116ee1cef2cSJed Brown
MatMultTransposeAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3)117d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N, Vec v1, Vec v2, Vec v3)
118d71ae5a4SJacob Faibussowitsch {
11954e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
120ee1cef2cSJed Brown
121ee1cef2cSJed Brown PetscFunctionBegin;
1229566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->lwork));
1239566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->lrestrict, v1, Na->lwork, INSERT_VALUES, SCATTER_REVERSE));
1249566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->lrestrict, v1, Na->lwork, INSERT_VALUES, SCATTER_REVERSE));
12511c5f74dSStefano Zampini if (v1 == v2) {
1269566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(Na->A, Na->lwork, Na->lwork, Na->rwork));
12711c5f74dSStefano Zampini } else if (v2 == v3) {
1289566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork));
1299566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, v2, Na->rwork, INSERT_VALUES, SCATTER_FORWARD));
1309566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, v2, Na->rwork, INSERT_VALUES, SCATTER_FORWARD));
1319566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(Na->A, Na->lwork, Na->rwork, Na->rwork));
1321e7f0bbdSJed Brown } else {
13311c5f74dSStefano Zampini if (!Na->rwork2) {
1349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->rwork, &Na->rwork2));
13511c5f74dSStefano Zampini } else {
1369566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Na->rwork2));
13711c5f74dSStefano Zampini }
1389566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, v2, Na->rwork2, INSERT_VALUES, SCATTER_FORWARD));
1399566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, v2, Na->rwork2, INSERT_VALUES, SCATTER_FORWARD));
1409566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(Na->A, Na->lwork, Na->rwork2, Na->rwork));
14111c5f74dSStefano Zampini }
1429566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Na->rprolong, Na->rwork, v3, INSERT_VALUES, SCATTER_REVERSE));
1439566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Na->rprolong, Na->rwork, v3, INSERT_VALUES, SCATTER_REVERSE));
1443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
145ee1cef2cSJed Brown }
146ee1cef2cSJed Brown
MatDestroy_SubMatrix(Mat N)147d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_SubMatrix(Mat N)
148d71ae5a4SJacob Faibussowitsch {
14954e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual *)N->data;
150ee1cef2cSJed Brown
151ee1cef2cSJed Brown PetscFunctionBegin;
1529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&Na->isrow));
1539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&Na->iscol));
1549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->lwork));
1559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rwork));
1569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->lwork2));
1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rwork2));
1589566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&Na->lrestrict));
1599566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&Na->rprolong));
1609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A));
1619566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data));
1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
163ee1cef2cSJed Brown }
164ee1cef2cSJed Brown
165ee1cef2cSJed Brown /*@
16611a5261eSBarry Smith MatCreateSubMatrixVirtual - Creates a virtual matrix `MATSUBMATRIX` that acts as a submatrix
167ee1cef2cSJed Brown
168c3339decSBarry Smith Collective
169ee1cef2cSJed Brown
170ee1cef2cSJed Brown Input Parameters:
171ee1cef2cSJed Brown + A - matrix that we will extract a submatrix of
172ee1cef2cSJed Brown . isrow - rows to be present in the submatrix
173ee1cef2cSJed Brown - iscol - columns to be present in the submatrix
174ee1cef2cSJed Brown
1752fe279fdSBarry Smith Output Parameter:
176ee1cef2cSJed Brown . newmat - new matrix
177ee1cef2cSJed Brown
178ee1cef2cSJed Brown Level: developer
179ee1cef2cSJed Brown
18011a5261eSBarry Smith Note:
18111a5261eSBarry Smith Most will use `MatCreateSubMatrix()` which provides a more efficient representation if it is available.
182ee1cef2cSJed Brown
1831cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSUBMATRIX`, `MATLOCALREF`, `MatCreateLocalRef()`, `MatCreateSubMatrix()`, `MatSubMatrixVirtualUpdate()`
184ee1cef2cSJed Brown @*/
MatCreateSubMatrixVirtual(Mat A,IS isrow,IS iscol,Mat * newmat)185d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrixVirtual(Mat A, IS isrow, IS iscol, Mat *newmat)
186d71ae5a4SJacob Faibussowitsch {
187ee1cef2cSJed Brown Vec left, right;
188ee1cef2cSJed Brown PetscInt m, n;
189ee1cef2cSJed Brown Mat N;
19054e05e6cSHong Zhang Mat_SubVirtual *Na;
191ee1cef2cSJed Brown
192ee1cef2cSJed Brown PetscFunctionBegin;
1930700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1940700a824SBarry Smith PetscValidHeaderSpecific(isrow, IS_CLASSID, 2);
1950700a824SBarry Smith PetscValidHeaderSpecific(iscol, IS_CLASSID, 3);
196*4f572ea9SToby Isaac PetscAssertPointer(newmat, 4);
197f4259b30SLisandro Dalcin *newmat = NULL;
198ee1cef2cSJed Brown
1999566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &N));
2009566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow, &m));
2019566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscol, &n));
2029566063dSJacob Faibussowitsch PetscCall(MatSetSizes(N, m, n, PETSC_DETERMINE, PETSC_DETERMINE));
2039566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)N, MATSUBMATRIX));
204ee1cef2cSJed Brown
2054dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&Na));
206ee1cef2cSJed Brown N->data = (void *)Na;
20711c5f74dSStefano Zampini
2089566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow));
2099566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol));
210ee1cef2cSJed Brown Na->isrow = isrow;
211ee1cef2cSJed Brown Na->iscol = iscol;
21211c5f74dSStefano Zampini
2139566063dSJacob Faibussowitsch PetscCall(PetscFree(N->defaultvectype));
2149566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(A->defaultvectype, &N->defaultvectype));
21511c5f74dSStefano Zampini /* Do not use MatConvert directly since MatShell has a duplicate operation which does not increase
21611c5f74dSStefano Zampini the reference count of the context. This is a problem if A is already of type MATSHELL */
2179566063dSJacob Faibussowitsch PetscCall(MatConvertFrom_Shell(A, MATSHELL, MAT_INITIAL_MATRIX, &Na->A));
218ee1cef2cSJed Brown
219ee1cef2cSJed Brown N->ops->destroy = MatDestroy_SubMatrix;
220ee1cef2cSJed Brown N->ops->mult = MatMult_SubMatrix;
221ee1cef2cSJed Brown N->ops->multadd = MatMultAdd_SubMatrix;
222ee1cef2cSJed Brown N->ops->multtranspose = MatMultTranspose_SubMatrix;
223ee1cef2cSJed Brown N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix;
224ee1cef2cSJed Brown N->ops->scale = MatScale_SubMatrix;
225ee1cef2cSJed Brown N->ops->diagonalscale = MatDiagonalScale_SubMatrix;
226bddd1ffdSAlp Dener N->ops->shift = MatShift_SubMatrix;
22711c5f74dSStefano Zampini N->ops->convert = MatConvert_Shell;
22811c5f74dSStefano Zampini N->ops->getdiagonal = MatGetDiagonal_SubMatrix;
229ee1cef2cSJed Brown
2309566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(N, A, A));
2319566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(N->rmap));
2329566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(N->cmap));
233ee1cef2cSJed Brown
2349566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &Na->rwork, &Na->lwork));
2359566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(N, &right, &left));
2369566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(Na->lwork, isrow, left, NULL, &Na->lrestrict));
2379566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(right, NULL, Na->rwork, iscol, &Na->rprolong));
2389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&left));
2399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&right));
2409566063dSJacob Faibussowitsch PetscCall(MatSetUp(N));
24126fbe8dcSKarl Rupp
24211c5f74dSStefano Zampini N->assembled = PETSC_TRUE;
243ee1cef2cSJed Brown *newmat = N;
2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
245ee1cef2cSJed Brown }
246ee1cef2cSJed Brown
24711a5261eSBarry Smith /*MC
24811a5261eSBarry Smith MATSUBMATRIX - "submatrix" - A matrix type that represents a virtual submatrix of a matrix
249ee1cef2cSJed Brown
25011a5261eSBarry Smith Level: advanced
25111a5261eSBarry Smith
25211a5261eSBarry Smith Developer Note:
2532ef1f0ffSBarry Smith The `MatType` is `MATSUBMATRIX` but the routines associated have `SubMatrixVirtual` in them, the `MatType` name should likely be changed to
2542ef1f0ffSBarry Smith `MATSUBMATRIXVIRTUAL`
25511a5261eSBarry Smith
2561cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateSubMatrixVirtual()`, `MatCreateSubMatrixVirtual()`, `MatCreateSubMatrix()`
25711a5261eSBarry Smith M*/
25811a5261eSBarry Smith
25911a5261eSBarry Smith /*@
26011a5261eSBarry Smith MatSubMatrixVirtualUpdate - Updates a `MATSUBMATRIX` virtual submatrix
26111a5261eSBarry Smith
262c3339decSBarry Smith Collective
263ee1cef2cSJed Brown
264ee1cef2cSJed Brown Input Parameters:
265ee1cef2cSJed Brown + N - submatrix to update
266ee1cef2cSJed Brown . A - full matrix in the submatrix
267ee1cef2cSJed Brown . isrow - rows in the update (same as the first time the submatrix was created)
268ee1cef2cSJed Brown - iscol - columns in the update (same as the first time the submatrix was created)
269ee1cef2cSJed Brown
270ee1cef2cSJed Brown Level: developer
271ee1cef2cSJed Brown
27211a5261eSBarry Smith Note:
27311a5261eSBarry Smith Most will use `MatCreateSubMatrix()` which provides a more efficient representation if it is available.
274ee1cef2cSJed Brown
2751cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATSUBMATRIX`, `MatCreateSubMatrixVirtual()`
276ee1cef2cSJed Brown @*/
MatSubMatrixVirtualUpdate(Mat N,Mat A,IS isrow,IS iscol)277d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSubMatrixVirtualUpdate(Mat N, Mat A, IS isrow, IS iscol)
278d71ae5a4SJacob Faibussowitsch {
279ace3abfcSBarry Smith PetscBool flg;
28054e05e6cSHong Zhang Mat_SubVirtual *Na;
281ee1cef2cSJed Brown
282ee1cef2cSJed Brown PetscFunctionBegin;
2830700a824SBarry Smith PetscValidHeaderSpecific(N, MAT_CLASSID, 1);
2840700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 2);
2850700a824SBarry Smith PetscValidHeaderSpecific(isrow, IS_CLASSID, 3);
2860700a824SBarry Smith PetscValidHeaderSpecific(iscol, IS_CLASSID, 4);
2879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)N, MATSUBMATRIX, &flg));
28828b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix has wrong type");
289ee1cef2cSJed Brown
29054e05e6cSHong Zhang Na = (Mat_SubVirtual *)N->data;
2919566063dSJacob Faibussowitsch PetscCall(ISEqual(isrow, Na->isrow, &flg));
29228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot update submatrix with different row indices");
2939566063dSJacob Faibussowitsch PetscCall(ISEqual(iscol, Na->iscol, &flg));
29428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot update submatrix with different column indices");
295ee1cef2cSJed Brown
2969566063dSJacob Faibussowitsch PetscCall(PetscFree(N->defaultvectype));
2979566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(A->defaultvectype, &N->defaultvectype));
2989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A));
29911c5f74dSStefano Zampini /* Do not use MatConvert directly since MatShell has a duplicate operation which does not increase
30011c5f74dSStefano Zampini the reference count of the context. This is a problem if A is already of type MATSHELL */
3019566063dSJacob Faibussowitsch PetscCall(MatConvertFrom_Shell(A, MATSHELL, MAT_INITIAL_MATRIX, &Na->A));
3023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
303ee1cef2cSJed Brown }
304