xref: /petsc/src/mat/impls/nest/matnest.c (revision b43b03e97385fe8ad94dfdfbf5044a665bcc78db)
1 
2 #include "matnestimpl.h" /*I   "petscmat.h"   I*/
3 
4 static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
5 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left);
6 
7 /* private functions */
8 #undef __FUNCT__
9 #define __FUNCT__ "MatNestGetSizes_Private"
10 static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
11 {
12   Mat_Nest       *bA = (Mat_Nest*)A->data;
13   PetscInt       i,j;
14   PetscErrorCode ierr;
15 
16   PetscFunctionBegin;
17   *m = *n = *M = *N = 0;
18   for (i=0; i<bA->nr; i++) {  /* rows */
19     PetscInt sm,sM;
20     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
21     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
22     *m += sm;
23     *M += sM;
24   }
25   for (j=0; j<bA->nc; j++) {  /* cols */
26     PetscInt sn,sN;
27     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
28     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
29     *n += sn;
30     *N += sN;
31   }
32   PetscFunctionReturn(0);
33 }
34 
35 /* operations */
36 #undef __FUNCT__
37 #define __FUNCT__ "MatMult_Nest"
38 static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
39 {
40   Mat_Nest       *bA = (Mat_Nest*)A->data;
41   Vec            *bx = bA->right,*by = bA->left;
42   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
43   PetscErrorCode ierr;
44 
45   PetscFunctionBegin;
46   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
47   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
48   for (i=0; i<nr; i++) {
49     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
50     for (j=0; j<nc; j++) {
51       if (!bA->m[i][j]) continue;
52       /* y[i] <- y[i] + A[i][j] * x[j] */
53       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
54     }
55   }
56   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
57   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
58   PetscFunctionReturn(0);
59 }
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "MatMultAdd_Nest"
63 static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
64 {
65   Mat_Nest       *bA = (Mat_Nest*)A->data;
66   Vec            *bx = bA->right,*bz = bA->left;
67   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
68   PetscErrorCode ierr;
69 
70   PetscFunctionBegin;
71   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
72   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
73   for (i=0; i<nr; i++) {
74     if (y != z) {
75       Vec by;
76       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
77       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
78       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
79     }
80     for (j=0; j<nc; j++) {
81       if (!bA->m[i][j]) continue;
82       /* y[i] <- y[i] + A[i][j] * x[j] */
83       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
84     }
85   }
86   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
87   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "MatMultTranspose_Nest"
93 static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
94 {
95   Mat_Nest       *bA = (Mat_Nest*)A->data;
96   Vec            *bx = bA->left,*by = bA->right;
97   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
98   PetscErrorCode ierr;
99 
100   PetscFunctionBegin;
101   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
102   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
103   for (j=0; j<nc; j++) {
104     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
105     for (i=0; i<nr; i++) {
106       if (!bA->m[j][i]) continue;
107       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
108       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
109     }
110   }
111   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
112   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "MatMultTransposeAdd_Nest"
118 static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
119 {
120   Mat_Nest       *bA = (Mat_Nest*)A->data;
121   Vec            *bx = bA->left,*bz = bA->right;
122   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
123   PetscErrorCode ierr;
124 
125   PetscFunctionBegin;
126   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
127   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
128   for (j=0; j<nc; j++) {
129     if (y != z) {
130       Vec by;
131       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
132       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
133       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
134     }
135     for (i=0; i<nr; i++) {
136       if (!bA->m[j][i]) continue;
137       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
138       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
139     }
140   }
141   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
142   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
143   PetscFunctionReturn(0);
144 }
145 
146 #undef __FUNCT__
147 #define __FUNCT__ "MatNestDestroyISList"
148 static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
149 {
150   PetscErrorCode ierr;
151   IS             *lst = *list;
152   PetscInt       i;
153 
154   PetscFunctionBegin;
155   if (!lst) PetscFunctionReturn(0);
156   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
157   ierr = PetscFree(lst);CHKERRQ(ierr);
158   *list = PETSC_NULL;
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "MatDestroy_Nest"
164 static PetscErrorCode MatDestroy_Nest(Mat A)
165 {
166   Mat_Nest       *vs = (Mat_Nest*)A->data;
167   PetscInt       i,j;
168   PetscErrorCode ierr;
169 
170   PetscFunctionBegin;
171   /* release the matrices and the place holders */
172   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
173   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
174   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
175   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
176 
177   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
178   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
179 
180   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
181 
182   /* release the matrices and the place holders */
183   if (vs->m) {
184     for (i=0; i<vs->nr; i++) {
185       for (j=0; j<vs->nc; j++) {
186         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
187       }
188       ierr = PetscFree( vs->m[i] );CHKERRQ(ierr);
189     }
190     ierr = PetscFree(vs->m);CHKERRQ(ierr);
191   }
192   ierr = PetscFree(A->data);CHKERRQ(ierr);
193 
194   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "",0);CHKERRQ(ierr);
195   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C",   "",0);CHKERRQ(ierr);
196   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
197   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
198   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetISs_C",      "",0);CHKERRQ(ierr);
199   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetLocalISs_C", "",0);CHKERRQ(ierr);
200   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
201   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
202   PetscFunctionReturn(0);
203 }
204 
205 #undef __FUNCT__
206 #define __FUNCT__ "MatAssemblyBegin_Nest"
207 static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
208 {
209   Mat_Nest       *vs = (Mat_Nest*)A->data;
210   PetscInt       i,j;
211   PetscErrorCode ierr;
212 
213   PetscFunctionBegin;
214   for (i=0; i<vs->nr; i++) {
215     for (j=0; j<vs->nc; j++) {
216       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
217     }
218   }
219   PetscFunctionReturn(0);
220 }
221 
222 #undef __FUNCT__
223 #define __FUNCT__ "MatAssemblyEnd_Nest"
224 static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
225 {
226   Mat_Nest       *vs = (Mat_Nest*)A->data;
227   PetscInt       i,j;
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   for (i=0; i<vs->nr; i++) {
232     for (j=0; j<vs->nc; j++) {
233       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
234     }
235   }
236   PetscFunctionReturn(0);
237 }
238 
239 #undef __FUNCT__
240 #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
241 static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
242 {
243   PetscErrorCode ierr;
244   Mat_Nest       *vs = (Mat_Nest*)A->data;
245   PetscInt       j;
246   Mat            sub;
247 
248   PetscFunctionBegin;
249   sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */
250   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
251   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
252   *B = sub;
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
258 static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
259 {
260   PetscErrorCode ierr;
261   Mat_Nest       *vs = (Mat_Nest*)A->data;
262   PetscInt       i;
263   Mat            sub;
264 
265   PetscFunctionBegin;
266   sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */
267   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
268   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
269   *B = sub;
270   PetscFunctionReturn(0);
271 }
272 
273 #undef __FUNCT__
274 #define __FUNCT__ "MatNestFindIS"
275 static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
276 {
277   PetscErrorCode ierr;
278   PetscInt       i;
279   PetscBool      flg;
280 
281   PetscFunctionBegin;
282   PetscValidPointer(list,3);
283   PetscValidHeaderSpecific(is,IS_CLASSID,4);
284   PetscValidIntPointer(found,5);
285   *found = -1;
286   for (i=0; i<n; i++) {
287     if (!list[i]) continue;
288     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
289     if (flg) {
290       *found = i;
291       PetscFunctionReturn(0);
292     }
293   }
294   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "MatNestGetRow"
300 /* Get a block row as a new MatNest */
301 static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
302 {
303   Mat_Nest       *vs = (Mat_Nest*)A->data;
304   char           keyname[256];
305   PetscErrorCode ierr;
306 
307   PetscFunctionBegin;
308   *B = PETSC_NULL;
309   ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr);
310   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
311   if (*B) PetscFunctionReturn(0);
312 
313   ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
314   (*B)->assembled = A->assembled;
315   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
316   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
317   PetscFunctionReturn(0);
318 }
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "MatNestFindSubMat"
322 static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
323 {
324   Mat_Nest       *vs = (Mat_Nest*)A->data;
325   PetscErrorCode ierr;
326   PetscInt       row,col;
327   PetscBool      same,isFullCol;
328 
329   PetscFunctionBegin;
330   /* Check if full column space. This is a hack */
331   isFullCol = PETSC_FALSE;
332   ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
333   if (same) {
334     PetscInt n,first,step,i,an,am,afirst,astep;
335     ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
336     ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
337     isFullCol = PETSC_TRUE;
338     for (i=0,an=0; i<vs->nc; i++) {
339       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
340       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
341       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
342       an += am;
343     }
344     if (an != n) isFullCol = PETSC_FALSE;
345   }
346 
347   if (isFullCol) {
348     PetscInt row;
349     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
350     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
351   } else {
352     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
353     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
354     *B = vs->m[row][col];
355   }
356   PetscFunctionReturn(0);
357 }
358 
359 #undef __FUNCT__
360 #define __FUNCT__ "MatGetSubMatrix_Nest"
361 static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
362 {
363   PetscErrorCode ierr;
364   Mat_Nest       *vs = (Mat_Nest*)A->data;
365   Mat            sub;
366 
367   PetscFunctionBegin;
368   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
369   switch (reuse) {
370   case MAT_INITIAL_MATRIX:
371     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
372     *B = sub;
373     break;
374   case MAT_REUSE_MATRIX:
375     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
376     break;
377   case MAT_IGNORE_MATRIX:       /* Nothing to do */
378     break;
379   }
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
385 PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
386 {
387   PetscErrorCode ierr;
388   Mat_Nest       *vs = (Mat_Nest*)A->data;
389   Mat            sub;
390 
391   PetscFunctionBegin;
392   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
393   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
394   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
395   *B = sub;
396   PetscFunctionReturn(0);
397 }
398 
399 #undef __FUNCT__
400 #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
401 static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
402 {
403   PetscErrorCode ierr;
404   Mat_Nest       *vs = (Mat_Nest*)A->data;
405   Mat            sub;
406 
407   PetscFunctionBegin;
408   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
409   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
410   if (sub) {
411     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
412     ierr = MatDestroy(B);CHKERRQ(ierr);
413   }
414   PetscFunctionReturn(0);
415 }
416 
417 #undef __FUNCT__
418 #define __FUNCT__ "MatGetDiagonal_Nest"
419 static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
420 {
421   Mat_Nest       *bA = (Mat_Nest*)A->data;
422   PetscInt       i;
423   PetscErrorCode ierr;
424 
425   PetscFunctionBegin;
426   for (i=0; i<bA->nr; i++) {
427     Vec bv;
428     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
429     if (bA->m[i][i]) {
430       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
431     } else {
432       ierr = VecSet(bv,1.0);CHKERRQ(ierr);
433     }
434     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
435   }
436   PetscFunctionReturn(0);
437 }
438 
439 #undef __FUNCT__
440 #define __FUNCT__ "MatDiagonalScale_Nest"
441 static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
442 {
443   Mat_Nest       *bA = (Mat_Nest*)A->data;
444   Vec            bl,*br;
445   PetscInt       i,j;
446   PetscErrorCode ierr;
447 
448   PetscFunctionBegin;
449   ierr = PetscMalloc(bA->nc*sizeof(Vec),&br);CHKERRQ(ierr);
450   for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
451   for (i=0; i<bA->nr; i++) {
452     ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
453     for (j=0; j<bA->nc; j++) {
454       if (bA->m[i][j]) {
455         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
456       }
457     }
458     ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
459   }
460   for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
461   ierr = PetscFree(br);CHKERRQ(ierr);
462   PetscFunctionReturn(0);
463 }
464 
465 #undef __FUNCT__
466 #define __FUNCT__ "MatScale_Nest"
467 static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
468 {
469   Mat_Nest       *bA = (Mat_Nest*)A->data;
470   PetscInt       i,j;
471   PetscErrorCode ierr;
472 
473   PetscFunctionBegin;
474   for (i=0; i<bA->nr; i++) {
475     for (j=0; j<bA->nc; j++) {
476       if (bA->m[i][j]) {
477         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
478       }
479     }
480   }
481   PetscFunctionReturn(0);
482 }
483 
484 #undef __FUNCT__
485 #define __FUNCT__ "MatShift_Nest"
486 static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
487 {
488   Mat_Nest       *bA = (Mat_Nest*)A->data;
489   PetscInt       i;
490   PetscErrorCode ierr;
491 
492   PetscFunctionBegin;
493   for (i=0; i<bA->nr; i++) {
494     if (!bA->m[i][i]) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support for shifting an empty diagonal block, insert a matrix in block (%D,%D)",i,i);
495     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
496   }
497   PetscFunctionReturn(0);
498 }
499 
500 #undef __FUNCT__
501 #define __FUNCT__ "MatGetVecs_Nest"
502 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
503 {
504   Mat_Nest       *bA = (Mat_Nest*)A->data;
505   Vec            *L,*R;
506   MPI_Comm       comm;
507   PetscInt       i,j;
508   PetscErrorCode ierr;
509 
510   PetscFunctionBegin;
511   comm = ((PetscObject)A)->comm;
512   if (right) {
513     /* allocate R */
514     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
515     /* Create the right vectors */
516     for (j=0; j<bA->nc; j++) {
517       for (i=0; i<bA->nr; i++) {
518         if (bA->m[i][j]) {
519           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
520           break;
521         }
522       }
523       if (i==bA->nr) {
524         /* have an empty column */
525         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
526       }
527     }
528     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
529     /* hand back control to the nest vector */
530     for (j=0; j<bA->nc; j++) {
531       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
532     }
533     ierr = PetscFree(R);CHKERRQ(ierr);
534   }
535 
536   if (left) {
537     /* allocate L */
538     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
539     /* Create the left vectors */
540     for (i=0; i<bA->nr; i++) {
541       for (j=0; j<bA->nc; j++) {
542         if (bA->m[i][j]) {
543           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
544           break;
545         }
546       }
547       if (j==bA->nc) {
548         /* have an empty row */
549         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
550       }
551     }
552 
553     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
554     for (i=0; i<bA->nr; i++) {
555       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
556     }
557 
558     ierr = PetscFree(L);CHKERRQ(ierr);
559   }
560   PetscFunctionReturn(0);
561 }
562 
563 #undef __FUNCT__
564 #define __FUNCT__ "MatView_Nest"
565 static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
566 {
567   Mat_Nest       *bA = (Mat_Nest*)A->data;
568   PetscBool      isascii;
569   PetscInt       i,j;
570   PetscErrorCode ierr;
571 
572   PetscFunctionBegin;
573   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
574   if (isascii) {
575 
576     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
577     PetscViewerASCIIPushTab( viewer );    /* push0 */
578     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
579 
580     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
581     for (i=0; i<bA->nr; i++) {
582       for (j=0; j<bA->nc; j++) {
583         const MatType type;
584         char name[256] = "",prefix[256] = "";
585         PetscInt NR,NC;
586         PetscBool isNest = PETSC_FALSE;
587 
588         if (!bA->m[i][j]) {
589           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
590           continue;
591         }
592         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
593         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
594         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof name,"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
595         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof prefix,"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
596         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
597 
598         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
599 
600         if (isNest) {
601           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
602           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
603           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
604         }
605       }
606     }
607     PetscViewerASCIIPopTab(viewer);    /* pop0 */
608   }
609   PetscFunctionReturn(0);
610 }
611 
612 #undef __FUNCT__
613 #define __FUNCT__ "MatZeroEntries_Nest"
614 static PetscErrorCode MatZeroEntries_Nest(Mat A)
615 {
616   Mat_Nest       *bA = (Mat_Nest*)A->data;
617   PetscInt       i,j;
618   PetscErrorCode ierr;
619 
620   PetscFunctionBegin;
621   for (i=0; i<bA->nr; i++) {
622     for (j=0; j<bA->nc; j++) {
623       if (!bA->m[i][j]) continue;
624       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
625     }
626   }
627   PetscFunctionReturn(0);
628 }
629 
630 #undef __FUNCT__
631 #define __FUNCT__ "MatDuplicate_Nest"
632 static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
633 {
634   Mat_Nest       *bA = (Mat_Nest*)A->data;
635   Mat            *b;
636   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
637   PetscErrorCode ierr;
638 
639   PetscFunctionBegin;
640   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
641   for (i=0; i<nr; i++) {
642     for (j=0; j<nc; j++) {
643       if (bA->m[i][j]) {
644         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
645       } else {
646         b[i*nc+j] = PETSC_NULL;
647       }
648     }
649   }
650   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
651   /* Give the new MatNest exclusive ownership */
652   for (i=0; i<nr*nc; i++) {
653     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
654   }
655   ierr = PetscFree(b);CHKERRQ(ierr);
656 
657   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
658   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
659   PetscFunctionReturn(0);
660 }
661 
662 /* nest api */
663 EXTERN_C_BEGIN
664 #undef __FUNCT__
665 #define __FUNCT__ "MatNestGetSubMat_Nest"
666 PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
667 {
668   Mat_Nest *bA = (Mat_Nest*)A->data;
669   PetscFunctionBegin;
670   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
671   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
672   *mat = bA->m[idxm][jdxm];
673   PetscFunctionReturn(0);
674 }
675 EXTERN_C_END
676 
677 #undef __FUNCT__
678 #define __FUNCT__ "MatNestGetSubMat"
679 /*@C
680  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
681 
682  Not collective
683 
684  Input Parameters:
685 +   A  - nest matrix
686 .   idxm - index of the matrix within the nest matrix
687 -   jdxm - index of the matrix within the nest matrix
688 
689  Output Parameter:
690 .   sub - matrix at index idxm,jdxm within the nest matrix
691 
692  Level: developer
693 
694 .seealso: MatNestGetSize(), MatNestGetSubMats()
695 @*/
696 PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
697 {
698   PetscErrorCode ierr;
699 
700   PetscFunctionBegin;
701   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
702   PetscFunctionReturn(0);
703 }
704 
705 EXTERN_C_BEGIN
706 #undef __FUNCT__
707 #define __FUNCT__ "MatNestSetSubMat_Nest"
708 PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
709 {
710   Mat_Nest *bA = (Mat_Nest*)A->data;
711   PetscInt m,n,M,N,mi,ni,Mi,Ni;
712   PetscErrorCode ierr;
713 
714   PetscFunctionBegin;
715   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
716   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
717   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
718   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
719   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
720   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
721   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
722   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
723   if (M != Mi || N != Ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni);
724   if (m != mi || n != ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni);
725   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
726   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
727   bA->m[idxm][jdxm] = mat;
728   PetscFunctionReturn(0);
729 }
730 EXTERN_C_END
731 
732 #undef __FUNCT__
733 #define __FUNCT__ "MatNestSetSubMat"
734 /*@C
735  MatNestSetSubMat - Set a single submatrix in the nest matrix.
736 
737  Logically collective on the submatrix communicator
738 
739  Input Parameters:
740 +   A  - nest matrix
741 .   idxm - index of the matrix within the nest matrix
742 .   jdxm - index of the matrix within the nest matrix
743 -   sub - matrix at index idxm,jdxm within the nest matrix
744 
745  Notes:
746  The new submatrix must have the same size and communicator as that block of the nest.
747 
748  This increments the reference count of the submatrix.
749 
750  Level: developer
751 
752 .seealso: MatNestSetSubMats(), MatNestGetSubMat()
753 @*/
754 PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
755 {
756   PetscErrorCode ierr;
757 
758   PetscFunctionBegin;
759   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
760   PetscFunctionReturn(0);
761 }
762 
763 EXTERN_C_BEGIN
764 #undef __FUNCT__
765 #define __FUNCT__ "MatNestGetSubMats_Nest"
766 PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
767 {
768   Mat_Nest *bA = (Mat_Nest*)A->data;
769   PetscFunctionBegin;
770   if (M)   { *M   = bA->nr; }
771   if (N)   { *N   = bA->nc; }
772   if (mat) { *mat = bA->m;  }
773   PetscFunctionReturn(0);
774 }
775 EXTERN_C_END
776 
777 #undef __FUNCT__
778 #define __FUNCT__ "MatNestGetSubMats"
779 /*@C
780  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
781 
782  Not collective
783 
784  Input Parameters:
785 .   A  - nest matrix
786 
787  Output Parameter:
788 +   M - number of rows in the nest matrix
789 .   N - number of cols in the nest matrix
790 -   mat - 2d array of matrices
791 
792  Notes:
793 
794  The user should not free the array mat.
795 
796  Level: developer
797 
798 .seealso: MatNestGetSize(), MatNestGetSubMat()
799 @*/
800 PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
801 {
802   PetscErrorCode ierr;
803 
804   PetscFunctionBegin;
805   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
806   PetscFunctionReturn(0);
807 }
808 
809 EXTERN_C_BEGIN
810 #undef __FUNCT__
811 #define __FUNCT__ "MatNestGetSize_Nest"
812 PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
813 {
814   Mat_Nest  *bA = (Mat_Nest*)A->data;
815 
816   PetscFunctionBegin;
817   if (M) { *M  = bA->nr; }
818   if (N) { *N  = bA->nc; }
819   PetscFunctionReturn(0);
820 }
821 EXTERN_C_END
822 
823 #undef __FUNCT__
824 #define __FUNCT__ "MatNestGetSize"
825 /*@C
826  MatNestGetSize - Returns the size of the nest matrix.
827 
828  Not collective
829 
830  Input Parameters:
831 .   A  - nest matrix
832 
833  Output Parameter:
834 +   M - number of rows in the nested mat
835 -   N - number of cols in the nested mat
836 
837  Notes:
838 
839  Level: developer
840 
841 .seealso: MatNestGetSubMat(), MatNestGetSubMats()
842 @*/
843 PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
844 {
845   PetscErrorCode ierr;
846 
847   PetscFunctionBegin;
848   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
849   PetscFunctionReturn(0);
850 }
851 
852 #undef __FUNCT__
853 #define __FUNCT__ "MatNestGetISs_Nest"
854 PETSC_EXTERN_C PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
855 {
856   Mat_Nest       *vs = (Mat_Nest*)A->data;
857   PetscInt i;
858 
859   PetscFunctionBegin;
860   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
861   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
862   PetscFunctionReturn(0);
863 }
864 
865 #undef __FUNCT__
866 #define __FUNCT__ "MatNestGetISs"
867 /*@C
868  MatNestGetISs - Returns the index sets partitioning the row and column spaces
869 
870  Not collective
871 
872  Input Parameters:
873 .   A  - nest matrix
874 
875  Output Parameter:
876 +   rows - array of row index sets
877 -   cols - array of column index sets
878 
879  Level: advanced
880 
881  Notes:
882  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
883 
884 .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
885 @*/
886 PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
887 {
888   PetscErrorCode ierr;
889 
890   PetscFunctionBegin;
891   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
892   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
893   PetscFunctionReturn(0);
894 }
895 
896 #undef __FUNCT__
897 #define __FUNCT__ "MatNestGetLocalISs_Nest"
898 PETSC_EXTERN_C PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
899 {
900   Mat_Nest       *vs = (Mat_Nest*)A->data;
901   PetscInt i;
902 
903   PetscFunctionBegin;
904   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
905   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
906   PetscFunctionReturn(0);
907 }
908 
909 #undef __FUNCT__
910 #define __FUNCT__ "MatNestGetLocalISs"
911 /*@C
912  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
913 
914  Not collective
915 
916  Input Parameters:
917 .   A  - nest matrix
918 
919  Output Parameter:
920 +   rows - array of row index sets (or PETSC_NULL to ignore)
921 -   cols - array of column index sets (or PETSC_NULL to ignore)
922 
923  Level: advanced
924 
925  Notes:
926  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
927 
928 .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
929 @*/
930 PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
931 {
932   PetscErrorCode ierr;
933 
934   PetscFunctionBegin;
935   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
936   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
937   PetscFunctionReturn(0);
938 }
939 
940 EXTERN_C_BEGIN
941 #undef __FUNCT__
942 #define __FUNCT__ "MatNestSetVecType_Nest"
943 PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
944 {
945   PetscErrorCode ierr;
946   PetscBool      flg;
947 
948   PetscFunctionBegin;
949   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
950   /* In reality, this only distinguishes VECNEST and "other" */
951   A->ops->getvecs       = flg ? MatGetVecs_Nest       : 0;
952  PetscFunctionReturn(0);
953 }
954 EXTERN_C_END
955 
956 #undef __FUNCT__
957 #define __FUNCT__ "MatNestSetVecType"
958 /*@C
959  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
960 
961  Not collective
962 
963  Input Parameters:
964 +  A  - nest matrix
965 -  vtype - type to use for creating vectors
966 
967  Notes:
968 
969  Level: developer
970 
971 .seealso: MatGetVecs()
972 @*/
973 PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
974 {
975   PetscErrorCode ierr;
976 
977   PetscFunctionBegin;
978   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
979   PetscFunctionReturn(0);
980 }
981 
982 EXTERN_C_BEGIN
983 #undef __FUNCT__
984 #define __FUNCT__ "MatNestSetSubMats_Nest"
985 PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
986 {
987   Mat_Nest       *s = (Mat_Nest*)A->data;
988   PetscInt       i,j,m,n,M,N;
989   PetscErrorCode ierr;
990 
991   PetscFunctionBegin;
992   s->nr = nr;
993   s->nc = nc;
994 
995   /* Create space for submatrices */
996   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
997   for (i=0; i<nr; i++) {
998     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
999   }
1000   for (i=0; i<nr; i++) {
1001     for (j=0; j<nc; j++) {
1002       s->m[i][j] = a[i*nc+j];
1003       if (a[i*nc+j]) {
1004         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1005       }
1006     }
1007   }
1008 
1009   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1010 
1011   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
1012   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
1013   for (i=0; i<nr; i++) s->row_len[i]=-1;
1014   for (j=0; j<nc; j++) s->col_len[j]=-1;
1015 
1016   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1017 
1018   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1019   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1020   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
1021   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1022   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1023   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
1024 
1025   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1026   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1027 
1028   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
1029   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
1030   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
1031   PetscFunctionReturn(0);
1032 }
1033 EXTERN_C_END
1034 
1035 #undef __FUNCT__
1036 #define __FUNCT__ "MatNestSetSubMats"
1037 /*@
1038    MatNestSetSubMats - Sets the nested submatrices
1039 
1040    Collective on Mat
1041 
1042    Input Parameter:
1043 +  N - nested matrix
1044 .  nr - number of nested row blocks
1045 .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1046 .  nc - number of nested column blocks
1047 .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1048 -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1049 
1050    Level: advanced
1051 
1052 .seealso: MatCreateNest(), MATNEST
1053 @*/
1054 PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1055 {
1056   PetscErrorCode ierr;
1057   PetscInt i;
1058 
1059   PetscFunctionBegin;
1060   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1061   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1062   if (nr && is_row) {
1063     PetscValidPointer(is_row,3);
1064     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1065   }
1066   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
1067   if (nc && is_col) {
1068     PetscValidPointer(is_col,5);
1069     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1070   }
1071   if (nr*nc) PetscValidPointer(a,6);
1072   ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr);
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 #undef __FUNCT__
1077 #define __FUNCT__ "MatNestCreateAggregateL2G_Private"
1078 static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog,ISLocalToGlobalMapping *ltogb)
1079 {
1080   PetscErrorCode ierr;
1081   PetscBool flg;
1082   PetscInt i,j,m,mi,*ix;
1083 
1084   PetscFunctionBegin;
1085   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
1086     if (islocal[i]) {
1087       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
1088       flg = PETSC_TRUE;       /* We found a non-trivial entry */
1089     } else {
1090       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
1091     }
1092     m += mi;
1093   }
1094   if (flg) {
1095     ierr = PetscMalloc(m*sizeof(*ix),&ix);CHKERRQ(ierr);
1096     for (i=0,n=0; i<n; i++) {
1097       ISLocalToGlobalMapping smap = PETSC_NULL;
1098       VecScatter scat;
1099       IS isreq;
1100       Vec lvec,gvec;
1101       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
1102       Mat sub;
1103 
1104       if (sizeof (*x) != sizeof(PetscScalar)) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support when scalars smaller than integers");
1105       if (colflg) {
1106         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1107       } else {
1108         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1109       }
1110       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,PETSC_NULL);CHKERRQ(ierr);}
1111       if (islocal[i]) {
1112         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
1113       } else {
1114         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
1115       }
1116       for (j=0; j<mi; j++) ix[m+j] = j;
1117       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
1118       /*
1119         Now we need to extract the monolithic global indices that correspond to the given split global indices.
1120         In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
1121         The approach here is ugly because it uses VecScatter to move indices.
1122        */
1123       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
1124       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
1125       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
1126       ierr = VecScatterCreate(gvec,isreq,lvec,PETSC_NULL,&scat);CHKERRQ(ierr);
1127       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
1128       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
1129       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
1130       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1131       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1132       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
1133       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
1134       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
1135       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
1136       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1137       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
1138       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
1139       m += mi;
1140     }
1141     ierr = ISLocalToGlobalMappingCreate(((PetscObject)A)->comm,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
1142     *ltogb = PETSC_NULL;
1143   } else {
1144     *ltog = PETSC_NULL;
1145     *ltogb = PETSC_NULL;
1146   }
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 
1151 /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1152 /*
1153   nprocessors = NP
1154   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
1155        proc 0: => (g_0,h_0,)
1156        proc 1: => (g_1,h_1,)
1157        ...
1158        proc nprocs-1: => (g_NP-1,h_NP-1,)
1159 
1160             proc 0:                      proc 1:                    proc nprocs-1:
1161     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
1162 
1163             proc 0:
1164     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
1165             proc 1:
1166     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
1167 
1168             proc NP-1:
1169     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
1170 */
1171 #undef __FUNCT__
1172 #define __FUNCT__ "MatSetUp_NestIS_Private"
1173 static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1174 {
1175   Mat_Nest       *vs = (Mat_Nest*)A->data;
1176   PetscInt       i,j,offset,n,nsum,bs;
1177   PetscErrorCode ierr;
1178   Mat            sub;
1179 
1180   PetscFunctionBegin;
1181   ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr);
1182   ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr);
1183   if (is_row) { /* valid IS is passed in */
1184     /* refs on is[] are incremeneted */
1185     for (i=0; i<vs->nr; i++) {
1186       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
1187       vs->isglobal.row[i] = is_row[i];
1188     }
1189   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
1190     nsum = 0;
1191     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
1192       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1193       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
1194       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
1195       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1196       nsum += n;
1197     }
1198     offset = 0;
1199 #if defined(PETSC_HAVE_MPI_EXSCAN)
1200     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1201 #else
1202     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sorry but this code requires MPI_EXSCAN that doesn't exist on your machine's version of MPI, install a MPI2 with PETSc to get this functionality");
1203 #endif
1204     for (i=0; i<vs->nr; i++) {
1205       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1206       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
1207       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1208       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1209       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
1210       offset += n;
1211     }
1212   }
1213 
1214   if (is_col) { /* valid IS is passed in */
1215     /* refs on is[] are incremeneted */
1216     for (j=0; j<vs->nc; j++) {
1217       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
1218       vs->isglobal.col[j] = is_col[j];
1219     }
1220   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
1221     offset = A->cmap->rstart;
1222     nsum = 0;
1223     for (j=0; j<vs->nc; j++) {
1224       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1225       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
1226       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
1227       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1228       nsum += n;
1229     }
1230     offset = 0;
1231 #if defined(PETSC_HAVE_MPI_EXSCAN)
1232     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1233 #else
1234     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sorry but this code requires MPI_EXSCAN that doesn't exist on your machine's version of MPI, install a MPI2 with PETSc to get this functionality");
1235 #endif
1236     for (j=0; j<vs->nc; j++) {
1237       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1238       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
1239       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1240       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1241       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
1242       offset += n;
1243     }
1244   }
1245 
1246   /* Set up the local ISs */
1247   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
1248   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
1249   for (i=0,offset=0; i<vs->nr; i++) {
1250     IS                     isloc;
1251     ISLocalToGlobalMapping rmap = PETSC_NULL;
1252     PetscInt               nlocal,bs;
1253     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1254     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);}
1255     if (rmap) {
1256       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1257       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1258       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1259       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1260     } else {
1261       nlocal = 0;
1262       isloc  = PETSC_NULL;
1263     }
1264     vs->islocal.row[i] = isloc;
1265     offset += nlocal;
1266   }
1267   for (i=0,offset=0; i<vs->nc; i++) {
1268     IS                     isloc;
1269     ISLocalToGlobalMapping cmap = PETSC_NULL;
1270     PetscInt               nlocal,bs;
1271     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1272     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);}
1273     if (cmap) {
1274       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1275       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1276       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1277       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1278     } else {
1279       nlocal = 0;
1280       isloc  = PETSC_NULL;
1281     }
1282     vs->islocal.col[i] = isloc;
1283     offset += nlocal;
1284   }
1285 
1286   /* Set up the aggregate ISLocalToGlobalMapping */
1287   {
1288     ISLocalToGlobalMapping rmap,rmapb,cmap,cmapb;
1289     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap,&rmapb);CHKERRQ(ierr);
1290     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap,&cmapb);CHKERRQ(ierr);
1291     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
1292     if (rmapb && cmapb) {ierr = MatSetLocalToGlobalMappingBlock(A,rmapb,cmapb);CHKERRQ(ierr);}
1293     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
1294     ierr = ISLocalToGlobalMappingDestroy(&rmapb);CHKERRQ(ierr);
1295     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
1296     ierr = ISLocalToGlobalMappingDestroy(&cmapb);CHKERRQ(ierr);
1297   }
1298 
1299 #if defined(PETSC_USE_DEBUG)
1300   for (i=0; i<vs->nr; i++) {
1301     for (j=0; j<vs->nc; j++) {
1302       PetscInt m,n,M,N,mi,ni,Mi,Ni;
1303       Mat B = vs->m[i][j];
1304       if (!B) continue;
1305       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
1306       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
1307       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
1308       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
1309       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
1310       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1311       if (M != Mi || N != Ni) SETERRQ6(((PetscObject)sub)->comm,PETSC_ERR_ARG_INCOMP,"Global sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",M,N,i,j,Mi,Ni);
1312       if (m != mi || n != ni) SETERRQ6(((PetscObject)sub)->comm,PETSC_ERR_ARG_INCOMP,"Local sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",m,n,i,j,mi,ni);
1313     }
1314   }
1315 #endif
1316 
1317   /* Set A->assembled if all non-null blocks are currently assembled */
1318   for (i=0; i<vs->nr; i++) {
1319     for (j=0; j<vs->nc; j++) {
1320       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1321     }
1322   }
1323   A->assembled = PETSC_TRUE;
1324   PetscFunctionReturn(0);
1325 }
1326 
1327 #undef __FUNCT__
1328 #define __FUNCT__ "MatCreateNest"
1329 /*@
1330    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1331 
1332    Collective on Mat
1333 
1334    Input Parameter:
1335 +  comm - Communicator for the new Mat
1336 .  nr - number of nested row blocks
1337 .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1338 .  nc - number of nested column blocks
1339 .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1340 -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1341 
1342    Output Parameter:
1343 .  B - new matrix
1344 
1345    Level: advanced
1346 
1347 .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1348 @*/
1349 PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1350 {
1351   Mat            A;
1352   PetscErrorCode ierr;
1353 
1354   PetscFunctionBegin;
1355   *B = 0;
1356   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1357   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1358   ierr = MatSetUp(A);CHKERRQ(ierr);
1359   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1360   *B = A;
1361   PetscFunctionReturn(0);
1362 }
1363 
1364 /*MC
1365   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1366 
1367   Level: intermediate
1368 
1369   Notes:
1370   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1371   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1372   It is usually used with DMComposite and DMCreateMatrix()
1373 
1374 .seealso: MatCreate(), MatType, MatCreateNest()
1375 M*/
1376 EXTERN_C_BEGIN
1377 #undef __FUNCT__
1378 #define __FUNCT__ "MatCreate_Nest"
1379 PetscErrorCode MatCreate_Nest(Mat A)
1380 {
1381   Mat_Nest       *s;
1382   PetscErrorCode ierr;
1383 
1384   PetscFunctionBegin;
1385   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1386   A->data         = (void*)s;
1387   s->nr = s->nc   = -1;
1388   s->m            = PETSC_NULL;
1389 
1390   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1391   A->ops->mult                  = MatMult_Nest;
1392   A->ops->multadd               = MatMultAdd_Nest;
1393   A->ops->multtranspose         = MatMultTranspose_Nest;
1394   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1395   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1396   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1397   A->ops->zeroentries           = MatZeroEntries_Nest;
1398   A->ops->duplicate             = MatDuplicate_Nest;
1399   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1400   A->ops->destroy               = MatDestroy_Nest;
1401   A->ops->view                  = MatView_Nest;
1402   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1403   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1404   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1405   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1406   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1407   A->ops->scale                 = MatScale_Nest;
1408   A->ops->shift                 = MatShift_Nest;
1409 
1410   A->spptr        = 0;
1411   A->same_nonzero = PETSC_FALSE;
1412   A->assembled    = PETSC_FALSE;
1413 
1414   /* expose Nest api's */
1415   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1416   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C",   "MatNestSetSubMat_Nest",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1417   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1418   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1419   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetISs_C",      "MatNestGetISs_Nest",      MatNestGetISs_Nest);CHKERRQ(ierr);
1420   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetLocalISs_C", "MatNestGetLocalISs_Nest", MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1421   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1422   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1423 
1424   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1425   PetscFunctionReturn(0);
1426 }
1427 EXTERN_C_END
1428