xref: /petsc/src/mat/impls/submat/submat.c (revision 54e05e6ce97fc93ccc7b11530501e7c5277b12eb)
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         left,right;       /* optional scaling */
7ee1cef2cSJed Brown   Vec         olwork,orwork;    /* work vectors outside the scatters, only touched by PreScale and only created if needed*/
8ee1cef2cSJed Brown   Vec         lwork,rwork;      /* work vectors inside the scatters */
9ee1cef2cSJed Brown   VecScatter  lrestrict,rprolong;
10ee1cef2cSJed Brown   Mat         A;
11ee1cef2cSJed Brown   PetscScalar scale;
12*54e05e6cSHong Zhang } Mat_SubVirtual;
13ee1cef2cSJed Brown 
14ee1cef2cSJed Brown static PetscErrorCode PreScaleLeft(Mat N,Vec x,Vec *xx)
15ee1cef2cSJed Brown {
16*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
17ee1cef2cSJed Brown   PetscErrorCode ierr;
18ee1cef2cSJed Brown 
19ee1cef2cSJed Brown   PetscFunctionBegin;
20ee1cef2cSJed Brown   if (!Na->left) {
21ee1cef2cSJed Brown     *xx = x;
22ee1cef2cSJed Brown   } else {
23ee1cef2cSJed Brown     if (!Na->olwork) {
24ee1cef2cSJed Brown       ierr = VecDuplicate(Na->left,&Na->olwork);CHKERRQ(ierr);
25ee1cef2cSJed Brown     }
261e7f0bbdSJed Brown     ierr = VecPointwiseMult(Na->olwork,x,Na->left);CHKERRQ(ierr);
27ee1cef2cSJed Brown     *xx  = Na->olwork;
28ee1cef2cSJed Brown   }
29ee1cef2cSJed Brown   PetscFunctionReturn(0);
30ee1cef2cSJed Brown }
31ee1cef2cSJed Brown 
32ee1cef2cSJed Brown static PetscErrorCode PreScaleRight(Mat N,Vec x,Vec *xx)
33ee1cef2cSJed Brown {
34*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
35ee1cef2cSJed Brown   PetscErrorCode ierr;
36ee1cef2cSJed Brown 
37ee1cef2cSJed Brown   PetscFunctionBegin;
38ee1cef2cSJed Brown   if (!Na->right) {
39ee1cef2cSJed Brown     *xx = x;
40ee1cef2cSJed Brown   } else {
41ee1cef2cSJed Brown     if (!Na->orwork) {
42ee1cef2cSJed Brown       ierr = VecDuplicate(Na->right,&Na->orwork);CHKERRQ(ierr);
43ee1cef2cSJed Brown     }
441e7f0bbdSJed Brown     ierr = VecPointwiseMult(Na->orwork,x,Na->right);CHKERRQ(ierr);
45ee1cef2cSJed Brown     *xx  = Na->orwork;
46ee1cef2cSJed Brown   }
47ee1cef2cSJed Brown   PetscFunctionReturn(0);
48ee1cef2cSJed Brown }
49ee1cef2cSJed Brown 
50ee1cef2cSJed Brown static PetscErrorCode PostScaleLeft(Mat N,Vec x)
51ee1cef2cSJed Brown {
52*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
53ee1cef2cSJed Brown   PetscErrorCode ierr;
54ee1cef2cSJed Brown 
55ee1cef2cSJed Brown   PetscFunctionBegin;
56ee1cef2cSJed Brown   if (Na->left) {
57ee1cef2cSJed Brown     ierr = VecPointwiseMult(x,x,Na->left);CHKERRQ(ierr);
58ee1cef2cSJed Brown   }
59ee1cef2cSJed Brown   PetscFunctionReturn(0);
60ee1cef2cSJed Brown }
61ee1cef2cSJed Brown 
62ee1cef2cSJed Brown static PetscErrorCode PostScaleRight(Mat N,Vec x)
63ee1cef2cSJed Brown {
64*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
65ee1cef2cSJed Brown   PetscErrorCode ierr;
66ee1cef2cSJed Brown 
67ee1cef2cSJed Brown   PetscFunctionBegin;
68ee1cef2cSJed Brown   if (Na->right) {
69ee1cef2cSJed Brown     ierr = VecPointwiseMult(x,x,Na->right);CHKERRQ(ierr);
70ee1cef2cSJed Brown   }
71ee1cef2cSJed Brown   PetscFunctionReturn(0);
72ee1cef2cSJed Brown }
73ee1cef2cSJed Brown 
74ee1cef2cSJed Brown static PetscErrorCode MatScale_SubMatrix(Mat N,PetscScalar scale)
75ee1cef2cSJed Brown {
76*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
77ee1cef2cSJed Brown 
78ee1cef2cSJed Brown   PetscFunctionBegin;
79ee1cef2cSJed Brown   Na->scale *= scale;
80ee1cef2cSJed Brown   PetscFunctionReturn(0);
81ee1cef2cSJed Brown }
82ee1cef2cSJed Brown 
83ee1cef2cSJed Brown static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N,Vec left,Vec right)
84ee1cef2cSJed Brown {
85*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
86ee1cef2cSJed Brown   PetscErrorCode ierr;
87ee1cef2cSJed Brown 
88ee1cef2cSJed Brown   PetscFunctionBegin;
89ee1cef2cSJed Brown   if (left) {
90ee1cef2cSJed Brown     if (!Na->left) {
91ee1cef2cSJed Brown       ierr = VecDuplicate(left,&Na->left);CHKERRQ(ierr);
92ee1cef2cSJed Brown       ierr = VecCopy(left,Na->left);CHKERRQ(ierr);
93ee1cef2cSJed Brown     } else {
94ee1cef2cSJed Brown       ierr = VecPointwiseMult(Na->left,left,Na->left);CHKERRQ(ierr);
95ee1cef2cSJed Brown     }
96ee1cef2cSJed Brown   }
97ee1cef2cSJed Brown   if (right) {
98ee1cef2cSJed Brown     if (!Na->right) {
99ee1cef2cSJed Brown       ierr = VecDuplicate(right,&Na->right);CHKERRQ(ierr);
100ee1cef2cSJed Brown       ierr = VecCopy(right,Na->right);CHKERRQ(ierr);
101ee1cef2cSJed Brown     } else {
102ee1cef2cSJed Brown       ierr = VecPointwiseMult(Na->right,right,Na->right);CHKERRQ(ierr);
103ee1cef2cSJed Brown     }
104ee1cef2cSJed Brown   }
105ee1cef2cSJed Brown   PetscFunctionReturn(0);
106ee1cef2cSJed Brown }
107ee1cef2cSJed Brown 
108ee1cef2cSJed Brown static PetscErrorCode MatMult_SubMatrix(Mat N,Vec x,Vec y)
109ee1cef2cSJed Brown {
110*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
111b7784de6SSatish Balay   Vec            xx  = 0;
112ee1cef2cSJed Brown   PetscErrorCode ierr;
113ee1cef2cSJed Brown 
114ee1cef2cSJed Brown   PetscFunctionBegin;
115ee1cef2cSJed Brown   ierr = PreScaleRight(N,x,&xx);CHKERRQ(ierr);
116ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr);
117ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119ee1cef2cSJed Brown   ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr);
120ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
121ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122ee1cef2cSJed Brown   ierr = PostScaleLeft(N,y);CHKERRQ(ierr);
123ee1cef2cSJed Brown   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
124ee1cef2cSJed Brown   PetscFunctionReturn(0);
125ee1cef2cSJed Brown }
126ee1cef2cSJed Brown 
127ee1cef2cSJed Brown static PetscErrorCode MatMultAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3)
128ee1cef2cSJed Brown {
129*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
130b7784de6SSatish Balay   Vec            xx  = 0;
131ee1cef2cSJed Brown   PetscErrorCode ierr;
132ee1cef2cSJed Brown 
133ee1cef2cSJed Brown   PetscFunctionBegin;
134ee1cef2cSJed Brown   ierr = PreScaleRight(N,v1,&xx);CHKERRQ(ierr);
135ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr);
136ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
137ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
138ee1cef2cSJed Brown   ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr);
1391e7f0bbdSJed Brown   if (v2 == v3) {
140d4a378daSJed Brown     if (Na->scale == (PetscScalar)1.0 && !Na->left) {
1411e7f0bbdSJed Brown       ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1421e7f0bbdSJed Brown       ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,v3,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1431e7f0bbdSJed Brown     } else {
1441e7f0bbdSJed Brown       if (!Na->olwork) {ierr = VecDuplicate(v3,&Na->olwork);CHKERRQ(ierr);}
1451e7f0bbdSJed Brown       ierr = VecScatterBegin(Na->lrestrict,Na->lwork,Na->olwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1461e7f0bbdSJed Brown       ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,Na->olwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1471e7f0bbdSJed Brown       ierr = PostScaleLeft(N,Na->olwork);CHKERRQ(ierr);
1481e7f0bbdSJed Brown       ierr = VecAXPY(v3,Na->scale,Na->olwork);CHKERRQ(ierr);
1491e7f0bbdSJed Brown     }
1501e7f0bbdSJed Brown   } else {
151ee1cef2cSJed Brown     ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
152ee1cef2cSJed Brown     ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
153ee1cef2cSJed Brown     ierr = PostScaleLeft(N,v3);CHKERRQ(ierr);
154ee1cef2cSJed Brown     ierr = VecAYPX(v3,Na->scale,v2);CHKERRQ(ierr);
1551e7f0bbdSJed Brown   }
156ee1cef2cSJed Brown   PetscFunctionReturn(0);
157ee1cef2cSJed Brown }
158ee1cef2cSJed Brown 
159ee1cef2cSJed Brown static PetscErrorCode MatMultTranspose_SubMatrix(Mat N,Vec x,Vec y)
160ee1cef2cSJed Brown {
161*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
162b7784de6SSatish Balay   Vec            xx  = 0;
163ee1cef2cSJed Brown   PetscErrorCode ierr;
164ee1cef2cSJed Brown 
165ee1cef2cSJed Brown   PetscFunctionBegin;
166ee1cef2cSJed Brown   ierr = PreScaleLeft(N,x,&xx);CHKERRQ(ierr);
167ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr);
168ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
169ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
170ee1cef2cSJed Brown   ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr);
1711e7f0bbdSJed Brown   ierr = VecScatterBegin(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1721e7f0bbdSJed Brown   ierr = VecScatterEnd  (Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
173ee1cef2cSJed Brown   ierr = PostScaleRight(N,y);CHKERRQ(ierr);
174ee1cef2cSJed Brown   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
175ee1cef2cSJed Brown   PetscFunctionReturn(0);
176ee1cef2cSJed Brown }
177ee1cef2cSJed Brown 
178ee1cef2cSJed Brown static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3)
179ee1cef2cSJed Brown {
180*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
181b7784de6SSatish Balay   Vec            xx  = 0;
182ee1cef2cSJed Brown   PetscErrorCode ierr;
183ee1cef2cSJed Brown 
184ee1cef2cSJed Brown   PetscFunctionBegin;
185ee1cef2cSJed Brown   ierr = PreScaleLeft(N,v1,&xx);CHKERRQ(ierr);
186ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr);
187ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
188ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
189ee1cef2cSJed Brown   ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr);
1901e7f0bbdSJed Brown   if (v2 == v3) {
191d4a378daSJed Brown     if (Na->scale == (PetscScalar)1.0 && !Na->right) {
1921e7f0bbdSJed Brown       ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1931e7f0bbdSJed Brown       ierr = VecScatterEnd  (Na->rprolong,Na->rwork,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1941e7f0bbdSJed Brown     } else {
1951e7f0bbdSJed Brown       if (!Na->orwork) {ierr = VecDuplicate(v3,&Na->orwork);CHKERRQ(ierr);}
1961e7f0bbdSJed Brown       ierr = VecScatterBegin(Na->rprolong,Na->rwork,Na->orwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1971e7f0bbdSJed Brown       ierr = VecScatterEnd  (Na->rprolong,Na->rwork,Na->orwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1981e7f0bbdSJed Brown       ierr = PostScaleRight(N,Na->orwork);CHKERRQ(ierr);
1991e7f0bbdSJed Brown       ierr = VecAXPY(v3,Na->scale,Na->orwork);CHKERRQ(ierr);
2001e7f0bbdSJed Brown     }
2011e7f0bbdSJed Brown   } else {
2021e7f0bbdSJed Brown     ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2031e7f0bbdSJed Brown     ierr = VecScatterEnd  (Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
204ee1cef2cSJed Brown     ierr = PostScaleRight(N,v3);CHKERRQ(ierr);
205ee1cef2cSJed Brown     ierr = VecAYPX(v3,Na->scale,v2);CHKERRQ(ierr);
2061e7f0bbdSJed Brown   }
207ee1cef2cSJed Brown   PetscFunctionReturn(0);
208ee1cef2cSJed Brown }
209ee1cef2cSJed Brown 
210ee1cef2cSJed Brown static PetscErrorCode MatDestroy_SubMatrix(Mat N)
211ee1cef2cSJed Brown {
212*54e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
213ee1cef2cSJed Brown   PetscErrorCode ierr;
214ee1cef2cSJed Brown 
215ee1cef2cSJed Brown   PetscFunctionBegin;
2166bf464f9SBarry Smith   ierr = ISDestroy(&Na->isrow);CHKERRQ(ierr);
2176bf464f9SBarry Smith   ierr = ISDestroy(&Na->iscol);CHKERRQ(ierr);
2186bf464f9SBarry Smith   ierr = VecDestroy(&Na->left);CHKERRQ(ierr);
2196bf464f9SBarry Smith   ierr = VecDestroy(&Na->right);CHKERRQ(ierr);
2206bf464f9SBarry Smith   ierr = VecDestroy(&Na->olwork);CHKERRQ(ierr);
2216bf464f9SBarry Smith   ierr = VecDestroy(&Na->orwork);CHKERRQ(ierr);
2226bf464f9SBarry Smith   ierr = VecDestroy(&Na->lwork);CHKERRQ(ierr);
2236bf464f9SBarry Smith   ierr = VecDestroy(&Na->rwork);CHKERRQ(ierr);
2246bf464f9SBarry Smith   ierr = VecScatterDestroy(&Na->lrestrict);CHKERRQ(ierr);
2256bf464f9SBarry Smith   ierr = VecScatterDestroy(&Na->rprolong);CHKERRQ(ierr);
2266bf464f9SBarry Smith   ierr = MatDestroy(&Na->A);CHKERRQ(ierr);
2276bf464f9SBarry Smith   ierr = PetscFree(N->data);CHKERRQ(ierr);
228ee1cef2cSJed Brown   PetscFunctionReturn(0);
229ee1cef2cSJed Brown }
230ee1cef2cSJed Brown 
231ee1cef2cSJed Brown /*@
232*54e05e6cSHong Zhang    MatCreateSubMatrixVirtual - Creates a virtual matrix that acts as a submatrix
233ee1cef2cSJed Brown 
234ee1cef2cSJed Brown    Collective on Mat
235ee1cef2cSJed Brown 
236ee1cef2cSJed Brown    Input Parameters:
237ee1cef2cSJed Brown +  A - matrix that we will extract a submatrix of
238ee1cef2cSJed Brown .  isrow - rows to be present in the submatrix
239ee1cef2cSJed Brown -  iscol - columns to be present in the submatrix
240ee1cef2cSJed Brown 
241ee1cef2cSJed Brown    Output Parameters:
242ee1cef2cSJed Brown .  newmat - new matrix
243ee1cef2cSJed Brown 
244ee1cef2cSJed Brown    Level: developer
245ee1cef2cSJed Brown 
246ee1cef2cSJed Brown    Notes:
2477dae84e0SHong Zhang    Most will use MatCreateSubMatrix which provides a more efficient representation if it is available.
248ee1cef2cSJed Brown 
249*54e05e6cSHong Zhang .seealso: MatCreateSubMatrix(), MatSubMatrixVirtualUpdate()
250ee1cef2cSJed Brown @*/
251*54e05e6cSHong Zhang PetscErrorCode MatCreateSubMatrixVirtual(Mat A,IS isrow,IS iscol,Mat *newmat)
252ee1cef2cSJed Brown {
253ee1cef2cSJed Brown   Vec            left,right;
254ee1cef2cSJed Brown   PetscInt       m,n;
255ee1cef2cSJed Brown   Mat            N;
256*54e05e6cSHong Zhang   Mat_SubVirtual *Na;
257ee1cef2cSJed Brown   PetscErrorCode ierr;
258ee1cef2cSJed Brown 
259ee1cef2cSJed Brown   PetscFunctionBegin;
2600700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
2610700a824SBarry Smith   PetscValidHeaderSpecific(isrow,IS_CLASSID,2);
2620700a824SBarry Smith   PetscValidHeaderSpecific(iscol,IS_CLASSID,3);
263ee1cef2cSJed Brown   PetscValidPointer(newmat,4);
264ee1cef2cSJed Brown   *newmat = 0;
265ee1cef2cSJed Brown 
266ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&N);CHKERRQ(ierr);
267ee1cef2cSJed Brown   ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr);
268ee1cef2cSJed Brown   ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
269ee1cef2cSJed Brown   ierr = MatSetSizes(N,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
270ee1cef2cSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)N,MATSUBMATRIX);CHKERRQ(ierr);
271ee1cef2cSJed Brown 
272b00a9115SJed Brown   ierr      = PetscNewLog(N,&Na);CHKERRQ(ierr);
273ee1cef2cSJed Brown   N->data   = (void*)Na;
274ee1cef2cSJed Brown   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
275ee1cef2cSJed Brown   ierr      = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
276ee1cef2cSJed Brown   ierr      = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
277ee1cef2cSJed Brown   Na->A     = A;
278ee1cef2cSJed Brown   Na->isrow = isrow;
279ee1cef2cSJed Brown   Na->iscol = iscol;
280ee1cef2cSJed Brown   Na->scale = 1.0;
281ee1cef2cSJed Brown 
282ee1cef2cSJed Brown   N->ops->destroy          = MatDestroy_SubMatrix;
283ee1cef2cSJed Brown   N->ops->mult             = MatMult_SubMatrix;
284ee1cef2cSJed Brown   N->ops->multadd          = MatMultAdd_SubMatrix;
285ee1cef2cSJed Brown   N->ops->multtranspose    = MatMultTranspose_SubMatrix;
286ee1cef2cSJed Brown   N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix;
287ee1cef2cSJed Brown   N->ops->scale            = MatScale_SubMatrix;
288ee1cef2cSJed Brown   N->ops->diagonalscale    = MatDiagonalScale_SubMatrix;
289ee1cef2cSJed Brown 
29033d57670SJed Brown   ierr = MatSetBlockSizesFromMats(N,A,A);CHKERRQ(ierr);
29126283091SBarry Smith   ierr = PetscLayoutSetUp(N->rmap);CHKERRQ(ierr);
29226283091SBarry Smith   ierr = PetscLayoutSetUp(N->cmap);CHKERRQ(ierr);
293ee1cef2cSJed Brown 
2942a7a6963SBarry Smith   ierr = MatCreateVecs(A,&Na->rwork,&Na->lwork);CHKERRQ(ierr);
295ce94432eSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)isrow),&left);CHKERRQ(ierr);
296ce94432eSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)iscol),&right);CHKERRQ(ierr);
297ee1cef2cSJed Brown   ierr = VecSetSizes(left,m,PETSC_DETERMINE);CHKERRQ(ierr);
298ee1cef2cSJed Brown   ierr = VecSetSizes(right,n,PETSC_DETERMINE);CHKERRQ(ierr);
299ee1cef2cSJed Brown   ierr = VecSetUp(left);CHKERRQ(ierr);
300ee1cef2cSJed Brown   ierr = VecSetUp(right);CHKERRQ(ierr);
3010298fd71SBarry Smith   ierr = VecScatterCreate(Na->lwork,isrow,left,NULL,&Na->lrestrict);CHKERRQ(ierr);
3020298fd71SBarry Smith   ierr = VecScatterCreate(right,NULL,Na->rwork,iscol,&Na->rprolong);CHKERRQ(ierr);
3036bf464f9SBarry Smith   ierr = VecDestroy(&left);CHKERRQ(ierr);
3046bf464f9SBarry Smith   ierr = VecDestroy(&right);CHKERRQ(ierr);
305ee1cef2cSJed Brown 
306be7c243fSJed Brown   N->assembled = PETSC_TRUE;
30726fbe8dcSKarl Rupp 
308be7c243fSJed Brown   ierr = MatSetUp(N);CHKERRQ(ierr);
30926fbe8dcSKarl Rupp 
310ee1cef2cSJed Brown   *newmat      = N;
311ee1cef2cSJed Brown   PetscFunctionReturn(0);
312ee1cef2cSJed Brown }
313ee1cef2cSJed Brown 
314ee1cef2cSJed Brown 
315ee1cef2cSJed Brown /*@
316*54e05e6cSHong Zhang    MatSubMatrixVirtualUpdate - Updates a submatrix
317ee1cef2cSJed Brown 
318ee1cef2cSJed Brown    Collective on Mat
319ee1cef2cSJed Brown 
320ee1cef2cSJed Brown    Input Parameters:
321ee1cef2cSJed Brown +  N - submatrix to update
322ee1cef2cSJed Brown .  A - full matrix in the submatrix
323ee1cef2cSJed Brown .  isrow - rows in the update (same as the first time the submatrix was created)
324ee1cef2cSJed Brown -  iscol - columns in the update (same as the first time the submatrix was created)
325ee1cef2cSJed Brown 
326ee1cef2cSJed Brown    Level: developer
327ee1cef2cSJed Brown 
328ee1cef2cSJed Brown    Notes:
3297dae84e0SHong Zhang    Most will use MatCreateSubMatrix which provides a more efficient representation if it is available.
330ee1cef2cSJed Brown 
331*54e05e6cSHong Zhang .seealso: MatCreateSubMatrixVirtual()
332ee1cef2cSJed Brown @*/
333*54e05e6cSHong Zhang PetscErrorCode  MatSubMatrixVirtualUpdate(Mat N,Mat A,IS isrow,IS iscol)
334ee1cef2cSJed Brown {
335ee1cef2cSJed Brown   PetscErrorCode ierr;
336ace3abfcSBarry Smith   PetscBool      flg;
337*54e05e6cSHong Zhang   Mat_SubVirtual *Na;
338ee1cef2cSJed Brown 
339ee1cef2cSJed Brown   PetscFunctionBegin;
3400700a824SBarry Smith   PetscValidHeaderSpecific(N,MAT_CLASSID,1);
3410700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,2);
3420700a824SBarry Smith   PetscValidHeaderSpecific(isrow,IS_CLASSID,3);
3430700a824SBarry Smith   PetscValidHeaderSpecific(iscol,IS_CLASSID,4);
344251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)N,MATSUBMATRIX,&flg);CHKERRQ(ierr);
345ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix has wrong type");
346ee1cef2cSJed Brown 
347*54e05e6cSHong Zhang   Na   = (Mat_SubVirtual*)N->data;
348ee1cef2cSJed Brown   ierr = ISEqual(isrow,Na->isrow,&flg);CHKERRQ(ierr);
349e32f2f54SBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different row indices");
350ee1cef2cSJed Brown   ierr = ISEqual(iscol,Na->iscol,&flg);CHKERRQ(ierr);
351e32f2f54SBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different column indices");
352ee1cef2cSJed Brown 
353ee1cef2cSJed Brown   ierr  = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
3546bf464f9SBarry Smith   ierr  = MatDestroy(&Na->A);CHKERRQ(ierr);
355ee1cef2cSJed Brown   Na->A = A;
356ee1cef2cSJed Brown 
357ee1cef2cSJed Brown   Na->scale = 1.0;
3586bf464f9SBarry Smith   ierr      = VecDestroy(&Na->left);CHKERRQ(ierr);
3596bf464f9SBarry Smith   ierr      = VecDestroy(&Na->right);CHKERRQ(ierr);
360ee1cef2cSJed Brown   PetscFunctionReturn(0);
361ee1cef2cSJed Brown }
362