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