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