xref: /petsc/src/mat/impls/submat/submat.c (revision ee1cef2c4f92b4b45b8aa3ceac59cb944062fab9)
1*ee1cef2cSJed Brown #define PETSCMAT_DLL
2*ee1cef2cSJed Brown 
3*ee1cef2cSJed Brown #include "private/matimpl.h"          /*I "petscmat.h" I*/
4*ee1cef2cSJed Brown 
5*ee1cef2cSJed Brown typedef struct {
6*ee1cef2cSJed Brown   IS isrow,iscol;               /* rows and columns in submatrix, only used to check consistency */
7*ee1cef2cSJed Brown   Vec left,right;               /* optional scaling */
8*ee1cef2cSJed Brown   Vec olwork,orwork;            /* work vectors outside the scatters, only touched by PreScale and only created if needed*/
9*ee1cef2cSJed Brown   Vec lwork,rwork;              /* work vectors inside the scatters */
10*ee1cef2cSJed Brown   VecScatter lrestrict,rprolong;
11*ee1cef2cSJed Brown   Mat A;
12*ee1cef2cSJed Brown   PetscScalar scale;
13*ee1cef2cSJed Brown } Mat_SubMatrix;
14*ee1cef2cSJed Brown 
15*ee1cef2cSJed Brown #undef __FUNCT__
16*ee1cef2cSJed Brown #define __FUNCT__ "PreScaleLeft"
17*ee1cef2cSJed Brown static PetscErrorCode PreScaleLeft(Mat N,Vec x,Vec *xx)
18*ee1cef2cSJed Brown {
19*ee1cef2cSJed Brown   Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data;
20*ee1cef2cSJed Brown   PetscErrorCode ierr;
21*ee1cef2cSJed Brown 
22*ee1cef2cSJed Brown   PetscFunctionBegin;
23*ee1cef2cSJed Brown   if (!Na->left) {
24*ee1cef2cSJed Brown     *xx = x;
25*ee1cef2cSJed Brown   } else {
26*ee1cef2cSJed Brown     if (!Na->olwork) {
27*ee1cef2cSJed Brown       ierr = VecDuplicate(Na->left,&Na->olwork);CHKERRQ(ierr);
28*ee1cef2cSJed Brown     }
29*ee1cef2cSJed Brown     ierr = VecPointwiseMult(Na->left,x,Na->olwork);CHKERRQ(ierr);
30*ee1cef2cSJed Brown     *xx = Na->olwork;
31*ee1cef2cSJed Brown   }
32*ee1cef2cSJed Brown   PetscFunctionReturn(0);
33*ee1cef2cSJed Brown }
34*ee1cef2cSJed Brown 
35*ee1cef2cSJed Brown #undef __FUNCT__
36*ee1cef2cSJed Brown #define __FUNCT__ "PreScaleRight"
37*ee1cef2cSJed Brown static PetscErrorCode PreScaleRight(Mat N,Vec x,Vec *xx)
38*ee1cef2cSJed Brown {
39*ee1cef2cSJed Brown   Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data;
40*ee1cef2cSJed Brown   PetscErrorCode ierr;
41*ee1cef2cSJed Brown 
42*ee1cef2cSJed Brown   PetscFunctionBegin;
43*ee1cef2cSJed Brown   if (!Na->right) {
44*ee1cef2cSJed Brown     *xx = x;
45*ee1cef2cSJed Brown   } else {
46*ee1cef2cSJed Brown     if (!Na->orwork) {
47*ee1cef2cSJed Brown       ierr = VecDuplicate(Na->right,&Na->orwork);CHKERRQ(ierr);
48*ee1cef2cSJed Brown     }
49*ee1cef2cSJed Brown     ierr = VecPointwiseMult(Na->right,x,Na->orwork);CHKERRQ(ierr);
50*ee1cef2cSJed Brown     *xx = Na->orwork;
51*ee1cef2cSJed Brown   }
52*ee1cef2cSJed Brown   PetscFunctionReturn(0);
53*ee1cef2cSJed Brown }
54*ee1cef2cSJed Brown 
55*ee1cef2cSJed Brown #undef __FUNCT__
56*ee1cef2cSJed Brown #define __FUNCT__ "PostScaleLeft"
57*ee1cef2cSJed Brown static PetscErrorCode PostScaleLeft(Mat N,Vec x)
58*ee1cef2cSJed Brown {
59*ee1cef2cSJed Brown   Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data;
60*ee1cef2cSJed Brown   PetscErrorCode ierr;
61*ee1cef2cSJed Brown 
62*ee1cef2cSJed Brown   PetscFunctionBegin;
63*ee1cef2cSJed Brown   if (Na->left) {
64*ee1cef2cSJed Brown     ierr = VecPointwiseMult(x,x,Na->left);CHKERRQ(ierr);
65*ee1cef2cSJed Brown   }
66*ee1cef2cSJed Brown   PetscFunctionReturn(0);
67*ee1cef2cSJed Brown }
68*ee1cef2cSJed Brown 
69*ee1cef2cSJed Brown #undef __FUNCT__
70*ee1cef2cSJed Brown #define __FUNCT__ "PostScaleRight"
71*ee1cef2cSJed Brown static PetscErrorCode PostScaleRight(Mat N,Vec x)
72*ee1cef2cSJed Brown {
73*ee1cef2cSJed Brown   Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data;
74*ee1cef2cSJed Brown   PetscErrorCode ierr;
75*ee1cef2cSJed Brown 
76*ee1cef2cSJed Brown   PetscFunctionBegin;
77*ee1cef2cSJed Brown   if (Na->right) {
78*ee1cef2cSJed Brown     ierr = VecPointwiseMult(x,x,Na->right);CHKERRQ(ierr);
79*ee1cef2cSJed Brown   }
80*ee1cef2cSJed Brown   PetscFunctionReturn(0);
81*ee1cef2cSJed Brown }
82*ee1cef2cSJed Brown 
83*ee1cef2cSJed Brown #undef __FUNCT__
84*ee1cef2cSJed Brown #define __FUNCT__ "MatScale_SubMatrix"
85*ee1cef2cSJed Brown static PetscErrorCode MatScale_SubMatrix(Mat N,PetscScalar scale)
86*ee1cef2cSJed Brown {
87*ee1cef2cSJed Brown   Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data;
88*ee1cef2cSJed Brown 
89*ee1cef2cSJed Brown   PetscFunctionBegin;
90*ee1cef2cSJed Brown   Na->scale *= scale;
91*ee1cef2cSJed Brown   PetscFunctionReturn(0);
92*ee1cef2cSJed Brown }
93*ee1cef2cSJed Brown 
94*ee1cef2cSJed Brown #undef __FUNCT__
95*ee1cef2cSJed Brown #define __FUNCT__ "MatDiagonalScale_SubMatrix"
96*ee1cef2cSJed Brown static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N,Vec left,Vec right)
97*ee1cef2cSJed Brown {
98*ee1cef2cSJed Brown   Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data;
99*ee1cef2cSJed Brown   PetscErrorCode ierr;
100*ee1cef2cSJed Brown 
101*ee1cef2cSJed Brown   PetscFunctionBegin;
102*ee1cef2cSJed Brown   if (left) {
103*ee1cef2cSJed Brown     if (!Na->left) {
104*ee1cef2cSJed Brown       ierr = VecDuplicate(left,&Na->left);CHKERRQ(ierr);
105*ee1cef2cSJed Brown       ierr = VecCopy(left,Na->left);CHKERRQ(ierr);
106*ee1cef2cSJed Brown     } else {
107*ee1cef2cSJed Brown       ierr = VecPointwiseMult(Na->left,left,Na->left);CHKERRQ(ierr);
108*ee1cef2cSJed Brown     }
109*ee1cef2cSJed Brown   }
110*ee1cef2cSJed Brown   if (right) {
111*ee1cef2cSJed Brown     if (!Na->right) {
112*ee1cef2cSJed Brown       ierr = VecDuplicate(right,&Na->right);CHKERRQ(ierr);
113*ee1cef2cSJed Brown       ierr = VecCopy(right,Na->right);CHKERRQ(ierr);
114*ee1cef2cSJed Brown     } else {
115*ee1cef2cSJed Brown       ierr = VecPointwiseMult(Na->right,right,Na->right);CHKERRQ(ierr);
116*ee1cef2cSJed Brown     }
117*ee1cef2cSJed Brown   }
118*ee1cef2cSJed Brown   PetscFunctionReturn(0);
119*ee1cef2cSJed Brown }
120*ee1cef2cSJed Brown 
121*ee1cef2cSJed Brown #undef __FUNCT__
122*ee1cef2cSJed Brown #define __FUNCT__ "MatMult_SubMatrix"
123*ee1cef2cSJed Brown static PetscErrorCode MatMult_SubMatrix(Mat N,Vec x,Vec y)
124*ee1cef2cSJed Brown {
125*ee1cef2cSJed Brown   Mat_SubMatrix  *Na = (Mat_SubMatrix*)N->data;
126*ee1cef2cSJed Brown   Vec             xx;
127*ee1cef2cSJed Brown   PetscErrorCode  ierr;
128*ee1cef2cSJed Brown 
129*ee1cef2cSJed Brown   PetscFunctionBegin;
130*ee1cef2cSJed Brown   ierr = PreScaleRight(N,x,&xx);CHKERRQ(ierr);
131*ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr);
132*ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
133*ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
134*ee1cef2cSJed Brown   ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr);
135*ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
136*ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
137*ee1cef2cSJed Brown   ierr = PostScaleLeft(N,y);CHKERRQ(ierr);
138*ee1cef2cSJed Brown   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
139*ee1cef2cSJed Brown   PetscFunctionReturn(0);
140*ee1cef2cSJed Brown }
141*ee1cef2cSJed Brown 
142*ee1cef2cSJed Brown #undef __FUNCT__
143*ee1cef2cSJed Brown #define __FUNCT__ "MatMultAdd_SubMatrix"
144*ee1cef2cSJed Brown static PetscErrorCode MatMultAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3)
145*ee1cef2cSJed Brown {
146*ee1cef2cSJed Brown   Mat_SubMatrix  *Na = (Mat_SubMatrix*)N->data;
147*ee1cef2cSJed Brown   Vec             xx;
148*ee1cef2cSJed Brown   PetscErrorCode  ierr;
149*ee1cef2cSJed Brown 
150*ee1cef2cSJed Brown   PetscFunctionBegin;
151*ee1cef2cSJed Brown   ierr = PreScaleRight(N,v1,&xx);CHKERRQ(ierr);
152*ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr);
153*ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
154*ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
155*ee1cef2cSJed Brown   ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr);
156*ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
157*ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
158*ee1cef2cSJed Brown   ierr = PostScaleLeft(N,v3);CHKERRQ(ierr);
159*ee1cef2cSJed Brown   ierr = VecAYPX(v3,Na->scale,v2);CHKERRQ(ierr);
160*ee1cef2cSJed Brown   PetscFunctionReturn(0);
161*ee1cef2cSJed Brown }
162*ee1cef2cSJed Brown 
163*ee1cef2cSJed Brown #undef __FUNCT__
164*ee1cef2cSJed Brown #define __FUNCT__ "MatMultTranspose_SubMatrix"
165*ee1cef2cSJed Brown static PetscErrorCode MatMultTranspose_SubMatrix(Mat N,Vec x,Vec y)
166*ee1cef2cSJed Brown {
167*ee1cef2cSJed Brown   Mat_SubMatrix  *Na = (Mat_SubMatrix*)N->data;
168*ee1cef2cSJed Brown   Vec             xx;
169*ee1cef2cSJed Brown   PetscErrorCode ierr;
170*ee1cef2cSJed Brown 
171*ee1cef2cSJed Brown   PetscFunctionBegin;
172*ee1cef2cSJed Brown   ierr = PreScaleLeft(N,x,&xx);CHKERRQ(ierr);
173*ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr);
174*ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
175*ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
176*ee1cef2cSJed Brown   ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr);
177*ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
178*ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
179*ee1cef2cSJed Brown   ierr = PostScaleRight(N,y);CHKERRQ(ierr);
180*ee1cef2cSJed Brown   ierr = VecScale(y,Na->scale);CHKERRQ(ierr);
181*ee1cef2cSJed Brown   PetscFunctionReturn(0);
182*ee1cef2cSJed Brown }
183*ee1cef2cSJed Brown 
184*ee1cef2cSJed Brown #undef __FUNCT__
185*ee1cef2cSJed Brown #define __FUNCT__ "MatMultTransposeAdd_SubMatrix"
186*ee1cef2cSJed Brown static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3)
187*ee1cef2cSJed Brown {
188*ee1cef2cSJed Brown   Mat_SubMatrix  *Na = (Mat_SubMatrix*)N->data;
189*ee1cef2cSJed Brown   Vec             xx;
190*ee1cef2cSJed Brown   PetscErrorCode ierr;
191*ee1cef2cSJed Brown 
192*ee1cef2cSJed Brown   PetscFunctionBegin;
193*ee1cef2cSJed Brown   ierr = PreScaleLeft(N,v1,&xx);CHKERRQ(ierr);
194*ee1cef2cSJed Brown   ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr);
195*ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
196*ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
197*ee1cef2cSJed Brown   ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr);
198*ee1cef2cSJed Brown   ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
199*ee1cef2cSJed Brown   ierr = VecScatterEnd  (Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
200*ee1cef2cSJed Brown   ierr = PostScaleRight(N,v3);CHKERRQ(ierr);
201*ee1cef2cSJed Brown   ierr = VecAYPX(v3,Na->scale,v2);CHKERRQ(ierr);
202*ee1cef2cSJed Brown   PetscFunctionReturn(0);
203*ee1cef2cSJed Brown }
204*ee1cef2cSJed Brown 
205*ee1cef2cSJed Brown #undef __FUNCT__
206*ee1cef2cSJed Brown #define __FUNCT__ "MatDestroy_SubMatrix"
207*ee1cef2cSJed Brown static PetscErrorCode MatDestroy_SubMatrix(Mat N)
208*ee1cef2cSJed Brown {
209*ee1cef2cSJed Brown   Mat_SubMatrix  *Na = (Mat_SubMatrix*)N->data;
210*ee1cef2cSJed Brown   PetscErrorCode ierr;
211*ee1cef2cSJed Brown 
212*ee1cef2cSJed Brown   PetscFunctionBegin;
213*ee1cef2cSJed Brown   ierr = ISDestroy(Na->isrow);CHKERRQ(ierr);
214*ee1cef2cSJed Brown   ierr = ISDestroy(Na->iscol);CHKERRQ(ierr);
215*ee1cef2cSJed Brown   if (Na->left) {ierr = VecDestroy(Na->left);CHKERRQ(ierr);}
216*ee1cef2cSJed Brown   if (Na->right) {ierr = VecDestroy(Na->right);CHKERRQ(ierr);}
217*ee1cef2cSJed Brown   if (Na->olwork) {ierr = VecDestroy(Na->olwork);CHKERRQ(ierr);}
218*ee1cef2cSJed Brown   if (Na->orwork) {ierr = VecDestroy(Na->orwork);CHKERRQ(ierr);}
219*ee1cef2cSJed Brown   ierr = VecDestroy(Na->lwork);CHKERRQ(ierr);
220*ee1cef2cSJed Brown   ierr = VecDestroy(Na->rwork);CHKERRQ(ierr);
221*ee1cef2cSJed Brown   ierr = VecScatterDestroy(Na->lrestrict);CHKERRQ(ierr);
222*ee1cef2cSJed Brown   ierr = VecScatterDestroy(Na->rprolong);CHKERRQ(ierr);
223*ee1cef2cSJed Brown   ierr = MatDestroy(Na->A);CHKERRQ(ierr);
224*ee1cef2cSJed Brown   ierr = PetscFree(Na);CHKERRQ(ierr);
225*ee1cef2cSJed Brown   PetscFunctionReturn(0);
226*ee1cef2cSJed Brown }
227*ee1cef2cSJed Brown 
228*ee1cef2cSJed Brown #undef __FUNCT__
229*ee1cef2cSJed Brown #define __FUNCT__ "MatCreateSubMatrix"
230*ee1cef2cSJed Brown /*@
231*ee1cef2cSJed Brown    MatCreateSubMatrix - Creates a composite matrix that acts as a submatrix
232*ee1cef2cSJed Brown 
233*ee1cef2cSJed Brown    Collective on Mat
234*ee1cef2cSJed Brown 
235*ee1cef2cSJed Brown    Input Parameters:
236*ee1cef2cSJed Brown +  A - matrix that we will extract a submatrix of
237*ee1cef2cSJed Brown .  isrow - rows to be present in the submatrix
238*ee1cef2cSJed Brown -  iscol - columns to be present in the submatrix
239*ee1cef2cSJed Brown 
240*ee1cef2cSJed Brown    Output Parameters:
241*ee1cef2cSJed Brown .  newmat - new matrix
242*ee1cef2cSJed Brown 
243*ee1cef2cSJed Brown    Level: developer
244*ee1cef2cSJed Brown 
245*ee1cef2cSJed Brown    Notes:
246*ee1cef2cSJed Brown    Most will use MatGetSubMatrix which provides a more efficient representation if it is available.
247*ee1cef2cSJed Brown 
248*ee1cef2cSJed Brown .seealso: MatGetSubMatrix(), MatSubMatrixUpdate()
249*ee1cef2cSJed Brown @*/
250*ee1cef2cSJed Brown PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat A,IS isrow,IS iscol,Mat *newmat)
251*ee1cef2cSJed Brown {
252*ee1cef2cSJed Brown   Vec            left,right;
253*ee1cef2cSJed Brown   PetscInt       m,n;
254*ee1cef2cSJed Brown   Mat            N;
255*ee1cef2cSJed Brown   Mat_SubMatrix *Na;
256*ee1cef2cSJed Brown   PetscErrorCode ierr;
257*ee1cef2cSJed Brown 
258*ee1cef2cSJed Brown   PetscFunctionBegin;
259*ee1cef2cSJed Brown   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
260*ee1cef2cSJed Brown   PetscValidHeaderSpecific(isrow,IS_COOKIE,2);
261*ee1cef2cSJed Brown   PetscValidHeaderSpecific(iscol,IS_COOKIE,3);
262*ee1cef2cSJed Brown   PetscValidPointer(newmat,4);
263*ee1cef2cSJed Brown   *newmat = 0;
264*ee1cef2cSJed Brown 
265*ee1cef2cSJed Brown   ierr = MatCreate(((PetscObject)A)->comm,&N);CHKERRQ(ierr);
266*ee1cef2cSJed Brown   ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr);
267*ee1cef2cSJed Brown   ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
268*ee1cef2cSJed Brown   ierr = MatSetSizes(N,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
269*ee1cef2cSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)N,MATSUBMATRIX);CHKERRQ(ierr);
270*ee1cef2cSJed Brown 
271*ee1cef2cSJed Brown   ierr = PetscNewLog(N,Mat_SubMatrix,&Na);CHKERRQ(ierr);
272*ee1cef2cSJed Brown   N->data   = (void*)Na;
273*ee1cef2cSJed Brown   ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
274*ee1cef2cSJed Brown   ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
275*ee1cef2cSJed Brown   ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
276*ee1cef2cSJed Brown   Na->A     = A;
277*ee1cef2cSJed Brown   Na->isrow = isrow;
278*ee1cef2cSJed Brown   Na->iscol = iscol;
279*ee1cef2cSJed Brown   Na->scale = 1.0;
280*ee1cef2cSJed Brown 
281*ee1cef2cSJed Brown   N->ops->destroy          = MatDestroy_SubMatrix;
282*ee1cef2cSJed Brown   N->ops->mult             = MatMult_SubMatrix;
283*ee1cef2cSJed Brown   N->ops->multadd          = MatMultAdd_SubMatrix;
284*ee1cef2cSJed Brown   N->ops->multtranspose    = MatMultTranspose_SubMatrix;
285*ee1cef2cSJed Brown   N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix;
286*ee1cef2cSJed Brown   N->ops->scale            = MatScale_SubMatrix;
287*ee1cef2cSJed Brown   N->ops->diagonalscale    = MatDiagonalScale_SubMatrix;
288*ee1cef2cSJed Brown 
289*ee1cef2cSJed Brown   N->assembled = PETSC_TRUE;
290*ee1cef2cSJed Brown 
291*ee1cef2cSJed Brown   ierr = PetscMapSetBlockSize(N->rmap,A->rmap->bs);CHKERRQ(ierr);
292*ee1cef2cSJed Brown   ierr = PetscMapSetBlockSize(N->cmap,A->cmap->bs);CHKERRQ(ierr);
293*ee1cef2cSJed Brown   ierr = PetscMapSetUp(N->rmap);CHKERRQ(ierr);
294*ee1cef2cSJed Brown   ierr = PetscMapSetUp(N->cmap);CHKERRQ(ierr);
295*ee1cef2cSJed Brown 
296*ee1cef2cSJed Brown   ierr = MatGetVecs(A,&Na->rwork,&Na->lwork);CHKERRQ(ierr);
297*ee1cef2cSJed Brown   ierr = VecCreate(((PetscObject)isrow)->comm,&left);CHKERRQ(ierr);
298*ee1cef2cSJed Brown   ierr = VecCreate(((PetscObject)iscol)->comm,&right);CHKERRQ(ierr);
299*ee1cef2cSJed Brown   ierr = VecSetSizes(left,m,PETSC_DETERMINE);CHKERRQ(ierr);
300*ee1cef2cSJed Brown   ierr = VecSetSizes(right,n,PETSC_DETERMINE);CHKERRQ(ierr);
301*ee1cef2cSJed Brown   ierr = VecSetUp(left);CHKERRQ(ierr);
302*ee1cef2cSJed Brown   ierr = VecSetUp(right);CHKERRQ(ierr);
303*ee1cef2cSJed Brown   ierr = VecScatterCreate(Na->lwork,isrow,left,PETSC_NULL,&Na->lrestrict);CHKERRQ(ierr);
304*ee1cef2cSJed Brown   ierr = VecScatterCreate(right,PETSC_NULL,Na->rwork,iscol,&Na->rprolong);CHKERRQ(ierr);
305*ee1cef2cSJed Brown   ierr = VecDestroy(left);CHKERRQ(ierr);
306*ee1cef2cSJed Brown   ierr = VecDestroy(right);CHKERRQ(ierr);
307*ee1cef2cSJed Brown 
308*ee1cef2cSJed Brown   *newmat = N;
309*ee1cef2cSJed Brown   PetscFunctionReturn(0);
310*ee1cef2cSJed Brown }
311*ee1cef2cSJed Brown 
312*ee1cef2cSJed Brown 
313*ee1cef2cSJed Brown #undef __FUNCT__
314*ee1cef2cSJed Brown #define __FUNCT__ "MatSubMatrixUpdate"
315*ee1cef2cSJed Brown /*@
316*ee1cef2cSJed Brown    MatSubMatrixUpdate - Updates a submatrix
317*ee1cef2cSJed Brown 
318*ee1cef2cSJed Brown    Collective on Mat
319*ee1cef2cSJed Brown 
320*ee1cef2cSJed Brown    Input Parameters:
321*ee1cef2cSJed Brown +  N - submatrix to update
322*ee1cef2cSJed Brown .  A - full matrix in the submatrix
323*ee1cef2cSJed Brown .  isrow - rows in the update (same as the first time the submatrix was created)
324*ee1cef2cSJed Brown -  iscol - columns in the update (same as the first time the submatrix was created)
325*ee1cef2cSJed Brown 
326*ee1cef2cSJed Brown    Level: developer
327*ee1cef2cSJed Brown 
328*ee1cef2cSJed Brown    Notes:
329*ee1cef2cSJed Brown    Most will use MatGetSubMatrix which provides a more efficient representation if it is available.
330*ee1cef2cSJed Brown 
331*ee1cef2cSJed Brown .seealso: MatGetSubMatrix(), MatCreateSubMatrix()
332*ee1cef2cSJed Brown @*/
333*ee1cef2cSJed Brown PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat N,Mat A,IS isrow,IS iscol)
334*ee1cef2cSJed Brown {
335*ee1cef2cSJed Brown   PetscErrorCode  ierr;
336*ee1cef2cSJed Brown   PetscTruth      flg;
337*ee1cef2cSJed Brown   Mat_SubMatrix  *Na;
338*ee1cef2cSJed Brown 
339*ee1cef2cSJed Brown   PetscFunctionBegin;
340*ee1cef2cSJed Brown   PetscValidHeaderSpecific(N,MAT_COOKIE,1);
341*ee1cef2cSJed Brown   PetscValidHeaderSpecific(A,MAT_COOKIE,2);
342*ee1cef2cSJed Brown   PetscValidHeaderSpecific(isrow,IS_COOKIE,3);
343*ee1cef2cSJed Brown   PetscValidHeaderSpecific(iscol,IS_COOKIE,4);
344*ee1cef2cSJed Brown   ierr = PetscTypeCompare((PetscObject)N,MATSUBMATRIX,&flg);CHKERRQ(ierr);
345*ee1cef2cSJed Brown   if (!flg) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix has wrong type");
346*ee1cef2cSJed Brown 
347*ee1cef2cSJed Brown   Na = (Mat_SubMatrix*)N->data;
348*ee1cef2cSJed Brown   ierr = ISEqual(isrow,Na->isrow,&flg);CHKERRQ(ierr);
349*ee1cef2cSJed Brown   if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different row indices");
350*ee1cef2cSJed Brown   ierr = ISEqual(iscol,Na->iscol,&flg);CHKERRQ(ierr);
351*ee1cef2cSJed Brown   if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different column indices");
352*ee1cef2cSJed Brown 
353*ee1cef2cSJed Brown   ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
354*ee1cef2cSJed Brown   ierr = MatDestroy(Na->A);CHKERRQ(ierr);
355*ee1cef2cSJed Brown   Na->A = A;
356*ee1cef2cSJed Brown 
357*ee1cef2cSJed Brown   Na->scale = 1.0;
358*ee1cef2cSJed Brown   if (Na->left) {ierr = VecDestroy(Na->left);CHKERRQ(ierr);}
359*ee1cef2cSJed Brown   if (Na->right) {ierr = VecDestroy(Na->right);CHKERRQ(ierr);}
360*ee1cef2cSJed Brown   PetscFunctionReturn(0);
361*ee1cef2cSJed Brown }
362