xref: /petsc/src/mat/impls/submat/submat.c (revision bddd1ffd9e17d5d62698fccc2fcbd8387dd31d53)
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 */
9*bddd1ffdSAlp Dener   Vec         dshift;
10ee1cef2cSJed Brown   VecScatter  lrestrict,rprolong;
11ee1cef2cSJed Brown   Mat         A;
12*bddd1ffdSAlp Dener   PetscScalar vscale, axpy_vscale;
13*bddd1ffdSAlp Dener   PetscScalar vshift, axpy_vshift;
1454e05e6cSHong Zhang } Mat_SubVirtual;
15ee1cef2cSJed Brown 
16ee1cef2cSJed Brown static PetscErrorCode PreScaleLeft(Mat N,Vec x,Vec *xx)
17ee1cef2cSJed Brown {
1854e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
19ee1cef2cSJed Brown   PetscErrorCode ierr;
20ee1cef2cSJed Brown 
21ee1cef2cSJed Brown   PetscFunctionBegin;
22ee1cef2cSJed Brown   if (!Na->left) {
23ee1cef2cSJed Brown     *xx = x;
24ee1cef2cSJed Brown   } else {
25ee1cef2cSJed Brown     if (!Na->olwork) {
26ee1cef2cSJed Brown       ierr = VecDuplicate(Na->left,&Na->olwork);CHKERRQ(ierr);
27ee1cef2cSJed Brown     }
281e7f0bbdSJed Brown     ierr = VecPointwiseMult(Na->olwork,x,Na->left);CHKERRQ(ierr);
29ee1cef2cSJed Brown     *xx  = Na->olwork;
30ee1cef2cSJed Brown   }
31ee1cef2cSJed Brown   PetscFunctionReturn(0);
32ee1cef2cSJed Brown }
33ee1cef2cSJed Brown 
34ee1cef2cSJed Brown static PetscErrorCode PreScaleRight(Mat N,Vec x,Vec *xx)
35ee1cef2cSJed Brown {
3654e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
37ee1cef2cSJed Brown   PetscErrorCode ierr;
38ee1cef2cSJed Brown 
39ee1cef2cSJed Brown   PetscFunctionBegin;
40ee1cef2cSJed Brown   if (!Na->right) {
41ee1cef2cSJed Brown     *xx = x;
42ee1cef2cSJed Brown   } else {
43ee1cef2cSJed Brown     if (!Na->orwork) {
44ee1cef2cSJed Brown       ierr = VecDuplicate(Na->right,&Na->orwork);CHKERRQ(ierr);
45ee1cef2cSJed Brown     }
461e7f0bbdSJed Brown     ierr = VecPointwiseMult(Na->orwork,x,Na->right);CHKERRQ(ierr);
47ee1cef2cSJed Brown     *xx  = Na->orwork;
48ee1cef2cSJed Brown   }
49ee1cef2cSJed Brown   PetscFunctionReturn(0);
50ee1cef2cSJed Brown }
51ee1cef2cSJed Brown 
52ee1cef2cSJed Brown static PetscErrorCode PostScaleLeft(Mat N,Vec x)
53ee1cef2cSJed Brown {
5454e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
55ee1cef2cSJed Brown   PetscErrorCode ierr;
56ee1cef2cSJed Brown 
57ee1cef2cSJed Brown   PetscFunctionBegin;
58ee1cef2cSJed Brown   if (Na->left) {
59ee1cef2cSJed Brown     ierr = VecPointwiseMult(x,x,Na->left);CHKERRQ(ierr);
60ee1cef2cSJed Brown   }
61ee1cef2cSJed Brown   PetscFunctionReturn(0);
62ee1cef2cSJed Brown }
63ee1cef2cSJed Brown 
64ee1cef2cSJed Brown static PetscErrorCode PostScaleRight(Mat N,Vec x)
65ee1cef2cSJed Brown {
6654e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
67ee1cef2cSJed Brown   PetscErrorCode ierr;
68ee1cef2cSJed Brown 
69ee1cef2cSJed Brown   PetscFunctionBegin;
70ee1cef2cSJed Brown   if (Na->right) {
71ee1cef2cSJed Brown     ierr = VecPointwiseMult(x,x,Na->right);CHKERRQ(ierr);
72ee1cef2cSJed Brown   }
73ee1cef2cSJed Brown   PetscFunctionReturn(0);
74ee1cef2cSJed Brown }
75ee1cef2cSJed Brown 
76*bddd1ffdSAlp Dener /*
77*bddd1ffdSAlp Dener          Y = vscale*Y + diag(dshift)*X + vshift*X
78*bddd1ffdSAlp Dener 
79*bddd1ffdSAlp Dener          On input Y already contains A*x
80*bddd1ffdSAlp Dener */
81*bddd1ffdSAlp Dener static PetscErrorCode MatSubmatShiftAndScale(Mat A,Vec X,Vec Y)
82ee1cef2cSJed Brown {
83*bddd1ffdSAlp Dener   Mat_SubVirtual *Na = (Mat_SubVirtual*)A->data;
84*bddd1ffdSAlp Dener   PetscErrorCode ierr;
85ee1cef2cSJed Brown 
86ee1cef2cSJed Brown   PetscFunctionBegin;
87*bddd1ffdSAlp Dener   if (Na->dshift) {          /* get arrays because there is no VecPointwiseMultAdd() */
88*bddd1ffdSAlp Dener     PetscInt          i,m;
89*bddd1ffdSAlp Dener     const PetscScalar *x,*d;
90*bddd1ffdSAlp Dener     PetscScalar       *y;
91*bddd1ffdSAlp Dener     ierr = VecGetLocalSize(X,&m);CHKERRQ(ierr);
92*bddd1ffdSAlp Dener     ierr = VecGetArrayRead(Na->dshift,&d);CHKERRQ(ierr);
93*bddd1ffdSAlp Dener     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
94*bddd1ffdSAlp Dener     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
95*bddd1ffdSAlp Dener     for (i=0; i<m; i++) y[i] = Na->vscale*y[i] + d[i]*x[i];
96*bddd1ffdSAlp Dener     ierr = VecRestoreArrayRead(Na->dshift,&d);CHKERRQ(ierr);
97*bddd1ffdSAlp Dener     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
98*bddd1ffdSAlp Dener     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
99*bddd1ffdSAlp Dener   } else {
100*bddd1ffdSAlp Dener     ierr = VecScale(Y,Na->vscale);CHKERRQ(ierr);
101*bddd1ffdSAlp Dener   }
102*bddd1ffdSAlp Dener   if (Na->vshift != 0.0) {ierr = VecAXPY(Y,Na->vshift,X);CHKERRQ(ierr);} /* if test is for non-square matrices */
103*bddd1ffdSAlp Dener   PetscFunctionReturn(0);
104*bddd1ffdSAlp Dener }
105*bddd1ffdSAlp Dener 
106*bddd1ffdSAlp Dener static PetscErrorCode MatScale_SubMatrix(Mat N,PetscScalar a)
107*bddd1ffdSAlp Dener {
108*bddd1ffdSAlp Dener   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
109*bddd1ffdSAlp Dener   PetscErrorCode ierr;
110*bddd1ffdSAlp Dener 
111*bddd1ffdSAlp Dener   PetscFunctionBegin;
112*bddd1ffdSAlp Dener   Na->vscale *= a;
113*bddd1ffdSAlp Dener   Na->vshift *= a;
114*bddd1ffdSAlp Dener   if (Na->dshift) {
115*bddd1ffdSAlp Dener     ierr = VecScale(Na->dshift,a);CHKERRQ(ierr);
116*bddd1ffdSAlp Dener   }
117*bddd1ffdSAlp Dener   Na->axpy_vscale *= a;
118*bddd1ffdSAlp Dener   PetscFunctionReturn(0);
119*bddd1ffdSAlp Dener }
120*bddd1ffdSAlp Dener 
121*bddd1ffdSAlp Dener static PetscErrorCode MatShift_SubMatrix(Mat N,PetscScalar a)
122*bddd1ffdSAlp Dener {
123*bddd1ffdSAlp Dener   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
124*bddd1ffdSAlp Dener   PetscErrorCode ierr;
125*bddd1ffdSAlp Dener 
126*bddd1ffdSAlp Dener   PetscFunctionBegin;
127*bddd1ffdSAlp Dener   if (Na->left || Na->right) {
128*bddd1ffdSAlp Dener     if (!Na->dshift) {
129*bddd1ffdSAlp Dener       ierr = VecDuplicate(Na->left ? Na->left : Na->right, &Na->dshift);CHKERRQ(ierr);
130*bddd1ffdSAlp Dener       ierr = VecSet(Na->dshift,a);CHKERRQ(ierr);
131*bddd1ffdSAlp Dener     } else {
132*bddd1ffdSAlp Dener       if (Na->left)  {ierr = VecPointwiseMult(Na->dshift,Na->dshift,Na->left);CHKERRQ(ierr);}
133*bddd1ffdSAlp Dener       if (Na->right) {ierr = VecPointwiseMult(Na->dshift,Na->dshift,Na->right);CHKERRQ(ierr);}
134*bddd1ffdSAlp Dener       ierr = VecShift(Na->dshift,a);CHKERRQ(ierr);
135*bddd1ffdSAlp Dener     }
136*bddd1ffdSAlp Dener     if (Na->left)  {ierr = VecPointwiseDivide(Na->dshift,Na->dshift,Na->left);CHKERRQ(ierr);}
137*bddd1ffdSAlp Dener     if (Na->right) {ierr = VecPointwiseDivide(Na->dshift,Na->dshift,Na->right);CHKERRQ(ierr);}
138*bddd1ffdSAlp Dener   } else Na->vshift += a;
139ee1cef2cSJed Brown   PetscFunctionReturn(0);
140ee1cef2cSJed Brown }
141ee1cef2cSJed Brown 
142ee1cef2cSJed Brown static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N,Vec left,Vec right)
143ee1cef2cSJed Brown {
14454e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
145ee1cef2cSJed Brown   PetscErrorCode ierr;
146ee1cef2cSJed Brown 
147ee1cef2cSJed Brown   PetscFunctionBegin;
148ee1cef2cSJed Brown   if (left) {
149ee1cef2cSJed Brown     if (!Na->left) {
150ee1cef2cSJed Brown       ierr = VecDuplicate(left,&Na->left);CHKERRQ(ierr);
151ee1cef2cSJed Brown       ierr = VecCopy(left,Na->left);CHKERRQ(ierr);
152ee1cef2cSJed Brown     } else {
153ee1cef2cSJed Brown       ierr = VecPointwiseMult(Na->left,left,Na->left);CHKERRQ(ierr);
154ee1cef2cSJed Brown     }
155ee1cef2cSJed Brown   }
156ee1cef2cSJed Brown   if (right) {
157ee1cef2cSJed Brown     if (!Na->right) {
158ee1cef2cSJed Brown       ierr = VecDuplicate(right,&Na->right);CHKERRQ(ierr);
159ee1cef2cSJed Brown       ierr = VecCopy(right,Na->right);CHKERRQ(ierr);
160ee1cef2cSJed Brown     } else {
161ee1cef2cSJed Brown       ierr = VecPointwiseMult(Na->right,right,Na->right);CHKERRQ(ierr);
162ee1cef2cSJed Brown     }
163ee1cef2cSJed Brown   }
164ee1cef2cSJed Brown   PetscFunctionReturn(0);
165ee1cef2cSJed Brown }
166ee1cef2cSJed Brown 
167ee1cef2cSJed Brown static PetscErrorCode MatMult_SubMatrix(Mat N,Vec x,Vec y)
168ee1cef2cSJed Brown {
16954e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
170b7784de6SSatish Balay   Vec            xx  = 0;
171ee1cef2cSJed Brown   PetscErrorCode ierr;
172ee1cef2cSJed Brown 
173ee1cef2cSJed Brown   PetscFunctionBegin;
174ee1cef2cSJed Brown   ierr = PreScaleRight(N,x,&xx);CHKERRQ(ierr);
175ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr);
176ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
177ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
178ee1cef2cSJed Brown   ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr);
179ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
180ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
181*bddd1ffdSAlp Dener   ierr = MatSubmatShiftAndScale(N,xx,y);CHKERRQ(ierr);
182ee1cef2cSJed Brown   ierr = PostScaleLeft(N,y);CHKERRQ(ierr);
183ee1cef2cSJed Brown   PetscFunctionReturn(0);
184ee1cef2cSJed Brown }
185ee1cef2cSJed Brown 
186ee1cef2cSJed Brown static PetscErrorCode MatMultAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3)
187ee1cef2cSJed Brown {
18854e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
189b7784de6SSatish Balay   Vec            xx  = 0;
190ee1cef2cSJed Brown   PetscErrorCode ierr;
191ee1cef2cSJed Brown 
192ee1cef2cSJed Brown   PetscFunctionBegin;
193ee1cef2cSJed Brown   ierr = PreScaleRight(N,v1,&xx);CHKERRQ(ierr);
194ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr);
195ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
196ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
197ee1cef2cSJed Brown   ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr);
1981e7f0bbdSJed Brown   if (v2 == v3) {
1991e7f0bbdSJed Brown     if (!Na->olwork) {ierr = VecDuplicate(v3,&Na->olwork);CHKERRQ(ierr);}
2001e7f0bbdSJed Brown     ierr = VecScatterBegin(Na->lrestrict,Na->lwork,Na->olwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2011e7f0bbdSJed Brown     ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,Na->olwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
202*bddd1ffdSAlp Dener     ierr = MatSubmatShiftAndScale(N,xx,Na->olwork);CHKERRQ(ierr);
2031e7f0bbdSJed Brown     ierr = PostScaleLeft(N,Na->olwork);CHKERRQ(ierr);
204*bddd1ffdSAlp Dener     ierr = VecAXPY(v3,1.0,Na->olwork);CHKERRQ(ierr);
2051e7f0bbdSJed Brown   } else {
206ee1cef2cSJed Brown     ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
207ee1cef2cSJed Brown     ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
208*bddd1ffdSAlp Dener     ierr = MatSubmatShiftAndScale(N,xx,v3);CHKERRQ(ierr);
209ee1cef2cSJed Brown     ierr = PostScaleLeft(N,v3);CHKERRQ(ierr);
210*bddd1ffdSAlp Dener     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
2111e7f0bbdSJed Brown   }
212ee1cef2cSJed Brown   PetscFunctionReturn(0);
213ee1cef2cSJed Brown }
214ee1cef2cSJed Brown 
215ee1cef2cSJed Brown static PetscErrorCode MatMultTranspose_SubMatrix(Mat N,Vec x,Vec y)
216ee1cef2cSJed Brown {
21754e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
218b7784de6SSatish Balay   Vec            xx  = 0;
219ee1cef2cSJed Brown   PetscErrorCode ierr;
220ee1cef2cSJed Brown 
221ee1cef2cSJed Brown   PetscFunctionBegin;
222ee1cef2cSJed Brown   ierr = PreScaleLeft(N,x,&xx);CHKERRQ(ierr);
223ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr);
224ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
225ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
226ee1cef2cSJed Brown   ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr);
2271e7f0bbdSJed Brown   ierr = VecScatterBegin(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2281e7f0bbdSJed Brown   ierr = VecScatterEnd  (Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
229*bddd1ffdSAlp Dener   ierr = MatSubmatShiftAndScale(N,xx,y);CHKERRQ(ierr);
230ee1cef2cSJed Brown   ierr = PostScaleRight(N,y);CHKERRQ(ierr);
231ee1cef2cSJed Brown   PetscFunctionReturn(0);
232ee1cef2cSJed Brown }
233ee1cef2cSJed Brown 
234ee1cef2cSJed Brown static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3)
235ee1cef2cSJed Brown {
23654e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
237b7784de6SSatish Balay   Vec            xx  = 0;
238ee1cef2cSJed Brown   PetscErrorCode ierr;
239ee1cef2cSJed Brown 
240ee1cef2cSJed Brown   PetscFunctionBegin;
241ee1cef2cSJed Brown   ierr = PreScaleLeft(N,v1,&xx);CHKERRQ(ierr);
242ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr);
243ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
244ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
245ee1cef2cSJed Brown   ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr);
2461e7f0bbdSJed Brown   if (v2 == v3) {
2471e7f0bbdSJed Brown     if (!Na->orwork) {ierr = VecDuplicate(v3,&Na->orwork);CHKERRQ(ierr);}
2481e7f0bbdSJed Brown     ierr = VecScatterBegin(Na->rprolong,Na->rwork,Na->orwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2491e7f0bbdSJed Brown     ierr = VecScatterEnd  (Na->rprolong,Na->rwork,Na->orwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
250*bddd1ffdSAlp Dener     ierr = MatSubmatShiftAndScale(N,xx,Na->orwork);CHKERRQ(ierr);
2511e7f0bbdSJed Brown     ierr = PostScaleRight(N,Na->orwork);CHKERRQ(ierr);
252*bddd1ffdSAlp Dener     ierr = VecAXPY(v3,1.0,Na->orwork);CHKERRQ(ierr);
2531e7f0bbdSJed Brown   } else {
2541e7f0bbdSJed Brown     ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2551e7f0bbdSJed Brown     ierr = VecScatterEnd  (Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
256*bddd1ffdSAlp Dener     ierr = MatSubmatShiftAndScale(N,xx,v3);CHKERRQ(ierr);
257ee1cef2cSJed Brown     ierr = PostScaleRight(N,v3);CHKERRQ(ierr);
258*bddd1ffdSAlp Dener     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
2591e7f0bbdSJed Brown   }
260ee1cef2cSJed Brown   PetscFunctionReturn(0);
261ee1cef2cSJed Brown }
262ee1cef2cSJed Brown 
263ee1cef2cSJed Brown static PetscErrorCode MatDestroy_SubMatrix(Mat N)
264ee1cef2cSJed Brown {
26554e05e6cSHong Zhang   Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data;
266ee1cef2cSJed Brown   PetscErrorCode ierr;
267ee1cef2cSJed Brown 
268ee1cef2cSJed Brown   PetscFunctionBegin;
2696bf464f9SBarry Smith   ierr = ISDestroy(&Na->isrow);CHKERRQ(ierr);
2706bf464f9SBarry Smith   ierr = ISDestroy(&Na->iscol);CHKERRQ(ierr);
2716bf464f9SBarry Smith   ierr = VecDestroy(&Na->left);CHKERRQ(ierr);
2726bf464f9SBarry Smith   ierr = VecDestroy(&Na->right);CHKERRQ(ierr);
2736bf464f9SBarry Smith   ierr = VecDestroy(&Na->olwork);CHKERRQ(ierr);
2746bf464f9SBarry Smith   ierr = VecDestroy(&Na->orwork);CHKERRQ(ierr);
2756bf464f9SBarry Smith   ierr = VecDestroy(&Na->lwork);CHKERRQ(ierr);
2766bf464f9SBarry Smith   ierr = VecDestroy(&Na->rwork);CHKERRQ(ierr);
277*bddd1ffdSAlp Dener   ierr = VecDestroy(&Na->dshift);CHKERRQ(ierr);
2786bf464f9SBarry Smith   ierr = VecScatterDestroy(&Na->lrestrict);CHKERRQ(ierr);
2796bf464f9SBarry Smith   ierr = VecScatterDestroy(&Na->rprolong);CHKERRQ(ierr);
2806bf464f9SBarry Smith   ierr = MatDestroy(&Na->A);CHKERRQ(ierr);
2816bf464f9SBarry Smith   ierr = PetscFree(N->data);CHKERRQ(ierr);
282ee1cef2cSJed Brown   PetscFunctionReturn(0);
283ee1cef2cSJed Brown }
284ee1cef2cSJed Brown 
285ee1cef2cSJed Brown /*@
28654e05e6cSHong Zhang    MatCreateSubMatrixVirtual - Creates a virtual matrix that acts as a submatrix
287ee1cef2cSJed Brown 
288ee1cef2cSJed Brown    Collective on Mat
289ee1cef2cSJed Brown 
290ee1cef2cSJed Brown    Input Parameters:
291ee1cef2cSJed Brown +  A - matrix that we will extract a submatrix of
292ee1cef2cSJed Brown .  isrow - rows to be present in the submatrix
293ee1cef2cSJed Brown -  iscol - columns to be present in the submatrix
294ee1cef2cSJed Brown 
295ee1cef2cSJed Brown    Output Parameters:
296ee1cef2cSJed Brown .  newmat - new matrix
297ee1cef2cSJed Brown 
298ee1cef2cSJed Brown    Level: developer
299ee1cef2cSJed Brown 
300ee1cef2cSJed Brown    Notes:
3017dae84e0SHong Zhang    Most will use MatCreateSubMatrix which provides a more efficient representation if it is available.
302ee1cef2cSJed Brown 
30354e05e6cSHong Zhang .seealso: MatCreateSubMatrix(), MatSubMatrixVirtualUpdate()
304ee1cef2cSJed Brown @*/
30554e05e6cSHong Zhang PetscErrorCode MatCreateSubMatrixVirtual(Mat A,IS isrow,IS iscol,Mat *newmat)
306ee1cef2cSJed Brown {
307ee1cef2cSJed Brown   Vec            left,right;
308ee1cef2cSJed Brown   PetscInt       m,n;
309ee1cef2cSJed Brown   Mat            N;
31054e05e6cSHong Zhang   Mat_SubVirtual *Na;
311ee1cef2cSJed Brown   PetscErrorCode ierr;
312ee1cef2cSJed Brown 
313ee1cef2cSJed Brown   PetscFunctionBegin;
3140700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
3150700a824SBarry Smith   PetscValidHeaderSpecific(isrow,IS_CLASSID,2);
3160700a824SBarry Smith   PetscValidHeaderSpecific(iscol,IS_CLASSID,3);
317ee1cef2cSJed Brown   PetscValidPointer(newmat,4);
318ee1cef2cSJed Brown   *newmat = 0;
319ee1cef2cSJed Brown 
320ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&N);CHKERRQ(ierr);
321ee1cef2cSJed Brown   ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr);
322ee1cef2cSJed Brown   ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
323ee1cef2cSJed Brown   ierr = MatSetSizes(N,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
324ee1cef2cSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)N,MATSUBMATRIX);CHKERRQ(ierr);
325ee1cef2cSJed Brown 
326b00a9115SJed Brown   ierr      = PetscNewLog(N,&Na);CHKERRQ(ierr);
327ee1cef2cSJed Brown   N->data   = (void*)Na;
328ee1cef2cSJed Brown   ierr      = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
329ee1cef2cSJed Brown   ierr      = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
330ee1cef2cSJed Brown   ierr      = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
331ee1cef2cSJed Brown   Na->A     = A;
332ee1cef2cSJed Brown   Na->isrow = isrow;
333ee1cef2cSJed Brown   Na->iscol = iscol;
334*bddd1ffdSAlp Dener   Na->vscale = 1.0;
335*bddd1ffdSAlp Dener   Na->vshift = 0.0;
336ee1cef2cSJed Brown 
337ee1cef2cSJed Brown   N->ops->destroy          = MatDestroy_SubMatrix;
338ee1cef2cSJed Brown   N->ops->mult             = MatMult_SubMatrix;
339ee1cef2cSJed Brown   N->ops->multadd          = MatMultAdd_SubMatrix;
340ee1cef2cSJed Brown   N->ops->multtranspose    = MatMultTranspose_SubMatrix;
341ee1cef2cSJed Brown   N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix;
342ee1cef2cSJed Brown   N->ops->scale            = MatScale_SubMatrix;
343ee1cef2cSJed Brown   N->ops->diagonalscale    = MatDiagonalScale_SubMatrix;
344*bddd1ffdSAlp Dener   N->ops->shift            = MatShift_SubMatrix;
345ee1cef2cSJed Brown 
34633d57670SJed Brown   ierr = MatSetBlockSizesFromMats(N,A,A);CHKERRQ(ierr);
34726283091SBarry Smith   ierr = PetscLayoutSetUp(N->rmap);CHKERRQ(ierr);
34826283091SBarry Smith   ierr = PetscLayoutSetUp(N->cmap);CHKERRQ(ierr);
349ee1cef2cSJed Brown 
3502a7a6963SBarry Smith   ierr = MatCreateVecs(A,&Na->rwork,&Na->lwork);CHKERRQ(ierr);
351ce94432eSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)isrow),&left);CHKERRQ(ierr);
352ce94432eSBarry Smith   ierr = VecCreate(PetscObjectComm((PetscObject)iscol),&right);CHKERRQ(ierr);
353ee1cef2cSJed Brown   ierr = VecSetSizes(left,m,PETSC_DETERMINE);CHKERRQ(ierr);
354ee1cef2cSJed Brown   ierr = VecSetSizes(right,n,PETSC_DETERMINE);CHKERRQ(ierr);
355ee1cef2cSJed Brown   ierr = VecSetUp(left);CHKERRQ(ierr);
356ee1cef2cSJed Brown   ierr = VecSetUp(right);CHKERRQ(ierr);
3570298fd71SBarry Smith   ierr = VecScatterCreate(Na->lwork,isrow,left,NULL,&Na->lrestrict);CHKERRQ(ierr);
3580298fd71SBarry Smith   ierr = VecScatterCreate(right,NULL,Na->rwork,iscol,&Na->rprolong);CHKERRQ(ierr);
3596bf464f9SBarry Smith   ierr = VecDestroy(&left);CHKERRQ(ierr);
3606bf464f9SBarry Smith   ierr = VecDestroy(&right);CHKERRQ(ierr);
361ee1cef2cSJed Brown 
362be7c243fSJed Brown   N->assembled = PETSC_TRUE;
36326fbe8dcSKarl Rupp 
364be7c243fSJed Brown   ierr = MatSetUp(N);CHKERRQ(ierr);
36526fbe8dcSKarl Rupp 
366ee1cef2cSJed Brown   *newmat      = N;
367ee1cef2cSJed Brown   PetscFunctionReturn(0);
368ee1cef2cSJed Brown }
369ee1cef2cSJed Brown 
370ee1cef2cSJed Brown 
371ee1cef2cSJed Brown /*@
37254e05e6cSHong Zhang    MatSubMatrixVirtualUpdate - Updates a submatrix
373ee1cef2cSJed Brown 
374ee1cef2cSJed Brown    Collective on Mat
375ee1cef2cSJed Brown 
376ee1cef2cSJed Brown    Input Parameters:
377ee1cef2cSJed Brown +  N - submatrix to update
378ee1cef2cSJed Brown .  A - full matrix in the submatrix
379ee1cef2cSJed Brown .  isrow - rows in the update (same as the first time the submatrix was created)
380ee1cef2cSJed Brown -  iscol - columns in the update (same as the first time the submatrix was created)
381ee1cef2cSJed Brown 
382ee1cef2cSJed Brown    Level: developer
383ee1cef2cSJed Brown 
384ee1cef2cSJed Brown    Notes:
3857dae84e0SHong Zhang    Most will use MatCreateSubMatrix which provides a more efficient representation if it is available.
386ee1cef2cSJed Brown 
38754e05e6cSHong Zhang .seealso: MatCreateSubMatrixVirtual()
388ee1cef2cSJed Brown @*/
38954e05e6cSHong Zhang PetscErrorCode  MatSubMatrixVirtualUpdate(Mat N,Mat A,IS isrow,IS iscol)
390ee1cef2cSJed Brown {
391ee1cef2cSJed Brown   PetscErrorCode ierr;
392ace3abfcSBarry Smith   PetscBool      flg;
39354e05e6cSHong Zhang   Mat_SubVirtual *Na;
394ee1cef2cSJed Brown 
395ee1cef2cSJed Brown   PetscFunctionBegin;
3960700a824SBarry Smith   PetscValidHeaderSpecific(N,MAT_CLASSID,1);
3970700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,2);
3980700a824SBarry Smith   PetscValidHeaderSpecific(isrow,IS_CLASSID,3);
3990700a824SBarry Smith   PetscValidHeaderSpecific(iscol,IS_CLASSID,4);
400251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)N,MATSUBMATRIX,&flg);CHKERRQ(ierr);
401ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix has wrong type");
402ee1cef2cSJed Brown 
40354e05e6cSHong Zhang   Na   = (Mat_SubVirtual*)N->data;
404ee1cef2cSJed Brown   ierr = ISEqual(isrow,Na->isrow,&flg);CHKERRQ(ierr);
405e32f2f54SBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different row indices");
406ee1cef2cSJed Brown   ierr = ISEqual(iscol,Na->iscol,&flg);CHKERRQ(ierr);
407e32f2f54SBarry Smith   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different column indices");
408ee1cef2cSJed Brown 
409ee1cef2cSJed Brown   ierr  = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
4106bf464f9SBarry Smith   ierr  = MatDestroy(&Na->A);CHKERRQ(ierr);
411ee1cef2cSJed Brown   Na->A = A;
412ee1cef2cSJed Brown 
413*bddd1ffdSAlp Dener   Na->vshift = 0.0;
414*bddd1ffdSAlp Dener   Na->vscale = 1.0;
4156bf464f9SBarry Smith   ierr       = VecDestroy(&Na->left);CHKERRQ(ierr);
4166bf464f9SBarry Smith   ierr       = VecDestroy(&Na->right);CHKERRQ(ierr);
417ee1cef2cSJed Brown   PetscFunctionReturn(0);
418ee1cef2cSJed Brown }
419