xref: /petsc/src/mat/impls/nest/matnest.c (revision 4e6ef68f07cf4e7ffee1b93dbc438d90283d533b)
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 = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
469   if (r) {
470     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
471   }
472   bl = NULL;
473   for (i=0; i<bA->nr; i++) {
474     if (l) {
475       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
476     }
477     for (j=0; j<bA->nc; j++) {
478       if (bA->m[i][j]) {
479         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
480       }
481     }
482     if (l) {
483       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
484     }
485   }
486   if (r) {
487     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
488   }
489   ierr = PetscFree(br);CHKERRQ(ierr);
490   PetscFunctionReturn(0);
491 }
492 
493 #undef __FUNCT__
494 #define __FUNCT__ "MatScale_Nest"
495 static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
496 {
497   Mat_Nest       *bA = (Mat_Nest*)A->data;
498   PetscInt       i,j;
499   PetscErrorCode ierr;
500 
501   PetscFunctionBegin;
502   for (i=0; i<bA->nr; i++) {
503     for (j=0; j<bA->nc; j++) {
504       if (bA->m[i][j]) {
505         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
506       }
507     }
508   }
509   PetscFunctionReturn(0);
510 }
511 
512 #undef __FUNCT__
513 #define __FUNCT__ "MatShift_Nest"
514 static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
515 {
516   Mat_Nest       *bA = (Mat_Nest*)A->data;
517   PetscInt       i;
518   PetscErrorCode ierr;
519 
520   PetscFunctionBegin;
521   for (i=0; i<bA->nr; i++) {
522     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);
523     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
524   }
525   PetscFunctionReturn(0);
526 }
527 
528 #undef __FUNCT__
529 #define __FUNCT__ "MatGetVecs_Nest"
530 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
531 {
532   Mat_Nest       *bA = (Mat_Nest*)A->data;
533   Vec            *L,*R;
534   MPI_Comm       comm;
535   PetscInt       i,j;
536   PetscErrorCode ierr;
537 
538   PetscFunctionBegin;
539   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
540   if (right) {
541     /* allocate R */
542     ierr = PetscMalloc(sizeof(Vec) * bA->nc, &R);CHKERRQ(ierr);
543     /* Create the right vectors */
544     for (j=0; j<bA->nc; j++) {
545       for (i=0; i<bA->nr; i++) {
546         if (bA->m[i][j]) {
547           ierr = MatGetVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
548           break;
549         }
550       }
551       if (i==bA->nr) {
552         /* have an empty column */
553         SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
554       }
555     }
556     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
557     /* hand back control to the nest vector */
558     for (j=0; j<bA->nc; j++) {
559       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
560     }
561     ierr = PetscFree(R);CHKERRQ(ierr);
562   }
563 
564   if (left) {
565     /* allocate L */
566     ierr = PetscMalloc(sizeof(Vec) * bA->nr, &L);CHKERRQ(ierr);
567     /* Create the left vectors */
568     for (i=0; i<bA->nr; i++) {
569       for (j=0; j<bA->nc; j++) {
570         if (bA->m[i][j]) {
571           ierr = MatGetVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
572           break;
573         }
574       }
575       if (j==bA->nc) {
576         /* have an empty row */
577         SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
578       }
579     }
580 
581     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
582     for (i=0; i<bA->nr; i++) {
583       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
584     }
585 
586     ierr = PetscFree(L);CHKERRQ(ierr);
587   }
588   PetscFunctionReturn(0);
589 }
590 
591 #undef __FUNCT__
592 #define __FUNCT__ "MatView_Nest"
593 static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
594 {
595   Mat_Nest       *bA = (Mat_Nest*)A->data;
596   PetscBool      isascii;
597   PetscInt       i,j;
598   PetscErrorCode ierr;
599 
600   PetscFunctionBegin;
601   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
602   if (isascii) {
603 
604     PetscViewerASCIIPrintf(viewer,"Matrix object: \n");
605     PetscViewerASCIIPushTab(viewer);    /* push0 */
606     PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
607 
608     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");
609     for (i=0; i<bA->nr; i++) {
610       for (j=0; j<bA->nc; j++) {
611         MatType   type;
612         char      name[256] = "",prefix[256] = "";
613         PetscInt  NR,NC;
614         PetscBool isNest = PETSC_FALSE;
615 
616         if (!bA->m[i][j]) {
617           PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);
618           continue;
619         }
620         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
621         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
622         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
623         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
624         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
625 
626         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
627 
628         if (isNest) {
629           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
630           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
631           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
632         }
633       }
634     }
635     PetscViewerASCIIPopTab(viewer);    /* pop0 */
636   }
637   PetscFunctionReturn(0);
638 }
639 
640 #undef __FUNCT__
641 #define __FUNCT__ "MatZeroEntries_Nest"
642 static PetscErrorCode MatZeroEntries_Nest(Mat A)
643 {
644   Mat_Nest       *bA = (Mat_Nest*)A->data;
645   PetscInt       i,j;
646   PetscErrorCode ierr;
647 
648   PetscFunctionBegin;
649   for (i=0; i<bA->nr; i++) {
650     for (j=0; j<bA->nc; j++) {
651       if (!bA->m[i][j]) continue;
652       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
653     }
654   }
655   PetscFunctionReturn(0);
656 }
657 
658 #undef __FUNCT__
659 #define __FUNCT__ "MatDuplicate_Nest"
660 static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
661 {
662   Mat_Nest       *bA = (Mat_Nest*)A->data;
663   Mat            *b;
664   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
665   PetscErrorCode ierr;
666 
667   PetscFunctionBegin;
668   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
669   for (i=0; i<nr; i++) {
670     for (j=0; j<nc; j++) {
671       if (bA->m[i][j]) {
672         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
673       } else {
674         b[i*nc+j] = NULL;
675       }
676     }
677   }
678   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
679   /* Give the new MatNest exclusive ownership */
680   for (i=0; i<nr*nc; i++) {
681     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
682   }
683   ierr = PetscFree(b);CHKERRQ(ierr);
684 
685   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
686   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
687   PetscFunctionReturn(0);
688 }
689 
690 /* nest api */
691 #undef __FUNCT__
692 #define __FUNCT__ "MatNestGetSubMat_Nest"
693 PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
694 {
695   Mat_Nest *bA = (Mat_Nest*)A->data;
696 
697   PetscFunctionBegin;
698   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
699   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
700   *mat = bA->m[idxm][jdxm];
701   PetscFunctionReturn(0);
702 }
703 
704 #undef __FUNCT__
705 #define __FUNCT__ "MatNestGetSubMat"
706 /*@
707  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
708 
709  Not collective
710 
711  Input Parameters:
712 +   A  - nest matrix
713 .   idxm - index of the matrix within the nest matrix
714 -   jdxm - index of the matrix within the nest matrix
715 
716  Output Parameter:
717 .   sub - matrix at index idxm,jdxm within the nest matrix
718 
719  Level: developer
720 
721 .seealso: MatNestGetSize(), MatNestGetSubMats()
722 @*/
723 PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
724 {
725   PetscErrorCode ierr;
726 
727   PetscFunctionBegin;
728   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
729   PetscFunctionReturn(0);
730 }
731 
732 #undef __FUNCT__
733 #define __FUNCT__ "MatNestSetSubMat_Nest"
734 PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
735 {
736   Mat_Nest       *bA = (Mat_Nest*)A->data;
737   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
738   PetscErrorCode ierr;
739 
740   PetscFunctionBegin;
741   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
742   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
743   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
744   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
745   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
746   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
747   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
748   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
749   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);
750   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);
751 
752   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
753   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
754   bA->m[idxm][jdxm] = mat;
755   PetscFunctionReturn(0);
756 }
757 
758 #undef __FUNCT__
759 #define __FUNCT__ "MatNestSetSubMat"
760 /*@
761  MatNestSetSubMat - Set a single submatrix in the nest matrix.
762 
763  Logically collective on the submatrix communicator
764 
765  Input Parameters:
766 +   A  - nest matrix
767 .   idxm - index of the matrix within the nest matrix
768 .   jdxm - index of the matrix within the nest matrix
769 -   sub - matrix at index idxm,jdxm within the nest matrix
770 
771  Notes:
772  The new submatrix must have the same size and communicator as that block of the nest.
773 
774  This increments the reference count of the submatrix.
775 
776  Level: developer
777 
778 .seealso: MatNestSetSubMats(), MatNestGetSubMat()
779 @*/
780 PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
781 {
782   PetscErrorCode ierr;
783 
784   PetscFunctionBegin;
785   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
786   PetscFunctionReturn(0);
787 }
788 
789 #undef __FUNCT__
790 #define __FUNCT__ "MatNestGetSubMats_Nest"
791 PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
792 {
793   Mat_Nest *bA = (Mat_Nest*)A->data;
794 
795   PetscFunctionBegin;
796   if (M)   *M   = bA->nr;
797   if (N)   *N   = bA->nc;
798   if (mat) *mat = bA->m;
799   PetscFunctionReturn(0);
800 }
801 
802 #undef __FUNCT__
803 #define __FUNCT__ "MatNestGetSubMats"
804 /*@C
805  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
806 
807  Not collective
808 
809  Input Parameters:
810 .   A  - nest matrix
811 
812  Output Parameter:
813 +   M - number of rows in the nest matrix
814 .   N - number of cols in the nest matrix
815 -   mat - 2d array of matrices
816 
817  Notes:
818 
819  The user should not free the array mat.
820 
821  Level: developer
822 
823 .seealso: MatNestGetSize(), MatNestGetSubMat()
824 @*/
825 PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
826 {
827   PetscErrorCode ierr;
828 
829   PetscFunctionBegin;
830   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
831   PetscFunctionReturn(0);
832 }
833 
834 #undef __FUNCT__
835 #define __FUNCT__ "MatNestGetSize_Nest"
836 PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
837 {
838   Mat_Nest *bA = (Mat_Nest*)A->data;
839 
840   PetscFunctionBegin;
841   if (M) *M = bA->nr;
842   if (N) *N = bA->nc;
843   PetscFunctionReturn(0);
844 }
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "MatNestGetSize"
848 /*@
849  MatNestGetSize - Returns the size of the nest matrix.
850 
851  Not collective
852 
853  Input Parameters:
854 .   A  - nest matrix
855 
856  Output Parameter:
857 +   M - number of rows in the nested mat
858 -   N - number of cols in the nested mat
859 
860  Notes:
861 
862  Level: developer
863 
864 .seealso: MatNestGetSubMat(), MatNestGetSubMats()
865 @*/
866 PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
867 {
868   PetscErrorCode ierr;
869 
870   PetscFunctionBegin;
871   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
872   PetscFunctionReturn(0);
873 }
874 
875 #undef __FUNCT__
876 #define __FUNCT__ "MatNestGetISs_Nest"
877 static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
878 {
879   Mat_Nest *vs = (Mat_Nest*)A->data;
880   PetscInt i;
881 
882   PetscFunctionBegin;
883   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
884   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
885   PetscFunctionReturn(0);
886 }
887 
888 #undef __FUNCT__
889 #define __FUNCT__ "MatNestGetISs"
890 /*@C
891  MatNestGetISs - Returns the index sets partitioning the row and column spaces
892 
893  Not collective
894 
895  Input Parameters:
896 .   A  - nest matrix
897 
898  Output Parameter:
899 +   rows - array of row index sets
900 -   cols - array of column index sets
901 
902  Level: advanced
903 
904  Notes:
905  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
906 
907 .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
908 @*/
909 PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
910 {
911   PetscErrorCode ierr;
912 
913   PetscFunctionBegin;
914   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
915   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
916   PetscFunctionReturn(0);
917 }
918 
919 #undef __FUNCT__
920 #define __FUNCT__ "MatNestGetLocalISs_Nest"
921 static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
922 {
923   Mat_Nest *vs = (Mat_Nest*)A->data;
924   PetscInt i;
925 
926   PetscFunctionBegin;
927   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
928   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
929   PetscFunctionReturn(0);
930 }
931 
932 #undef __FUNCT__
933 #define __FUNCT__ "MatNestGetLocalISs"
934 /*@C
935  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
936 
937  Not collective
938 
939  Input Parameters:
940 .   A  - nest matrix
941 
942  Output Parameter:
943 +   rows - array of row index sets (or NULL to ignore)
944 -   cols - array of column index sets (or NULL to ignore)
945 
946  Level: advanced
947 
948  Notes:
949  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
950 
951 .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
952 @*/
953 PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
954 {
955   PetscErrorCode ierr;
956 
957   PetscFunctionBegin;
958   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
959   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
960   PetscFunctionReturn(0);
961 }
962 
963 #undef __FUNCT__
964 #define __FUNCT__ "MatNestSetVecType_Nest"
965 PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
966 {
967   PetscErrorCode ierr;
968   PetscBool      flg;
969 
970   PetscFunctionBegin;
971   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
972   /* In reality, this only distinguishes VECNEST and "other" */
973   if (flg) A->ops->getvecs = MatGetVecs_Nest;
974   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
975   PetscFunctionReturn(0);
976 }
977 
978 #undef __FUNCT__
979 #define __FUNCT__ "MatNestSetVecType"
980 /*@C
981  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
982 
983  Not collective
984 
985  Input Parameters:
986 +  A  - nest matrix
987 -  vtype - type to use for creating vectors
988 
989  Notes:
990 
991  Level: developer
992 
993 .seealso: MatGetVecs()
994 @*/
995 PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
996 {
997   PetscErrorCode ierr;
998 
999   PetscFunctionBegin;
1000   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1001   PetscFunctionReturn(0);
1002 }
1003 
1004 #undef __FUNCT__
1005 #define __FUNCT__ "MatNestSetSubMats_Nest"
1006 PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1007 {
1008   Mat_Nest       *s = (Mat_Nest*)A->data;
1009   PetscInt       i,j,m,n,M,N;
1010   PetscErrorCode ierr;
1011 
1012   PetscFunctionBegin;
1013   s->nr = nr;
1014   s->nc = nc;
1015 
1016   /* Create space for submatrices */
1017   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
1018   for (i=0; i<nr; i++) {
1019     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
1020   }
1021   for (i=0; i<nr; i++) {
1022     for (j=0; j<nc; j++) {
1023       s->m[i][j] = a[i*nc+j];
1024       if (a[i*nc+j]) {
1025         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1026       }
1027     }
1028   }
1029 
1030   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1031 
1032   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
1033   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
1034   for (i=0; i<nr; i++) s->row_len[i]=-1;
1035   for (j=0; j<nc; j++) s->col_len[j]=-1;
1036 
1037   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1038 
1039   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1040   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1041   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1042   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1043 
1044   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1045   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1046 
1047   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 #undef __FUNCT__
1052 #define __FUNCT__ "MatNestSetSubMats"
1053 /*@
1054    MatNestSetSubMats - Sets the nested submatrices
1055 
1056    Collective on Mat
1057 
1058    Input Parameter:
1059 +  N - nested matrix
1060 .  nr - number of nested row blocks
1061 .  is_row - index sets for each nested row block, or NULL to make contiguous
1062 .  nc - number of nested column blocks
1063 .  is_col - index sets for each nested column block, or NULL to make contiguous
1064 -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1065 
1066    Level: advanced
1067 
1068 .seealso: MatCreateNest(), MATNEST
1069 @*/
1070 PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1071 {
1072   PetscErrorCode ierr;
1073   PetscInt       i;
1074 
1075   PetscFunctionBegin;
1076   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1077   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1078   if (nr && is_row) {
1079     PetscValidPointer(is_row,3);
1080     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1081   }
1082   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
1083   if (nc && is_col) {
1084     PetscValidPointer(is_col,5);
1085     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1086   }
1087   if (nr*nc) PetscValidPointer(a,6);
1088   ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr);
1089   PetscFunctionReturn(0);
1090 }
1091 
1092 #undef __FUNCT__
1093 #define __FUNCT__ "MatNestCreateAggregateL2G_Private"
1094 static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog,ISLocalToGlobalMapping *ltogb)
1095 {
1096   PetscErrorCode ierr;
1097   PetscBool      flg;
1098   PetscInt       i,j,m,mi,*ix;
1099 
1100   PetscFunctionBegin;
1101   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
1102     if (islocal[i]) {
1103       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
1104       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
1105     } else {
1106       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
1107     }
1108     m += mi;
1109   }
1110   if (flg) {
1111     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
1112     for (i=0,n=0; i<n; i++) {
1113       ISLocalToGlobalMapping smap = NULL;
1114       VecScatter             scat;
1115       IS                     isreq;
1116       Vec                    lvec,gvec;
1117       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
1118       Mat sub;
1119 
1120       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
1121       if (colflg) {
1122         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1123       } else {
1124         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1125       }
1126       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);}
1127       if (islocal[i]) {
1128         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
1129       } else {
1130         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
1131       }
1132       for (j=0; j<mi; j++) ix[m+j] = j;
1133       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
1134       /*
1135         Now we need to extract the monolithic global indices that correspond to the given split global indices.
1136         In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
1137         The approach here is ugly because it uses VecScatter to move indices.
1138        */
1139       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
1140       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
1141       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
1142       ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr);
1143       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
1144       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
1145       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
1146       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1147       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1148       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
1149       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
1150       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
1151       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
1152       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1153       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
1154       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
1155       m   += mi;
1156     }
1157     ierr   = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
1158     *ltogb = NULL;
1159   } else {
1160     *ltog  = NULL;
1161     *ltogb = NULL;
1162   }
1163   PetscFunctionReturn(0);
1164 }
1165 
1166 
1167 /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1168 /*
1169   nprocessors = NP
1170   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1171        proc 0: => (g_0,h_0,)
1172        proc 1: => (g_1,h_1,)
1173        ...
1174        proc nprocs-1: => (g_NP-1,h_NP-1,)
1175 
1176             proc 0:                      proc 1:                    proc nprocs-1:
1177     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1178 
1179             proc 0:
1180     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1181             proc 1:
1182     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1183 
1184             proc NP-1:
1185     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1186 */
1187 #undef __FUNCT__
1188 #define __FUNCT__ "MatSetUp_NestIS_Private"
1189 static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1190 {
1191   Mat_Nest       *vs = (Mat_Nest*)A->data;
1192   PetscInt       i,j,offset,n,nsum,bs;
1193   PetscErrorCode ierr;
1194   Mat            sub = NULL;
1195 
1196   PetscFunctionBegin;
1197   ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr);
1198   ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr);
1199   if (is_row) { /* valid IS is passed in */
1200     /* refs on is[] are incremeneted */
1201     for (i=0; i<vs->nr; i++) {
1202       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
1203 
1204       vs->isglobal.row[i] = is_row[i];
1205     }
1206   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
1207     nsum = 0;
1208     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
1209       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1210       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
1211       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1212       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1213       nsum += n;
1214     }
1215     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1216     offset -= nsum;
1217     for (i=0; i<vs->nr; i++) {
1218       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1219       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1220       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1221       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1222       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
1223       offset += n;
1224     }
1225   }
1226 
1227   if (is_col) { /* valid IS is passed in */
1228     /* refs on is[] are incremeneted */
1229     for (j=0; j<vs->nc; j++) {
1230       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
1231 
1232       vs->isglobal.col[j] = is_col[j];
1233     }
1234   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
1235     offset = A->cmap->rstart;
1236     nsum   = 0;
1237     for (j=0; j<vs->nc; j++) {
1238       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1239       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
1240       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1241       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1242       nsum += n;
1243     }
1244     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1245     offset -= nsum;
1246     for (j=0; j<vs->nc; j++) {
1247       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1248       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1249       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1250       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1251       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
1252       offset += n;
1253     }
1254   }
1255 
1256   /* Set up the local ISs */
1257   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1258   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1259   for (i=0,offset=0; i<vs->nr; i++) {
1260     IS                     isloc;
1261     ISLocalToGlobalMapping rmap = NULL;
1262     PetscInt               nlocal,bs;
1263     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1264     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1265     if (rmap) {
1266       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1267       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1268       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1269       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1270     } else {
1271       nlocal = 0;
1272       isloc  = NULL;
1273     }
1274     vs->islocal.row[i] = isloc;
1275     offset            += nlocal;
1276   }
1277   for (i=0,offset=0; i<vs->nc; i++) {
1278     IS                     isloc;
1279     ISLocalToGlobalMapping cmap = NULL;
1280     PetscInt               nlocal,bs;
1281     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1282     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1283     if (cmap) {
1284       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1285       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1286       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1287       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1288     } else {
1289       nlocal = 0;
1290       isloc  = NULL;
1291     }
1292     vs->islocal.col[i] = isloc;
1293     offset            += nlocal;
1294   }
1295 
1296   /* Set up the aggregate ISLocalToGlobalMapping */
1297   {
1298     ISLocalToGlobalMapping rmap,rmapb,cmap,cmapb;
1299     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap,&rmapb);CHKERRQ(ierr);
1300     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap,&cmapb);CHKERRQ(ierr);
1301     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
1302     if (rmapb && cmapb) {ierr = MatSetLocalToGlobalMappingBlock(A,rmapb,cmapb);CHKERRQ(ierr);}
1303     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
1304     ierr = ISLocalToGlobalMappingDestroy(&rmapb);CHKERRQ(ierr);
1305     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
1306     ierr = ISLocalToGlobalMappingDestroy(&cmapb);CHKERRQ(ierr);
1307   }
1308 
1309 #if defined(PETSC_USE_DEBUG)
1310   for (i=0; i<vs->nr; i++) {
1311     for (j=0; j<vs->nc; j++) {
1312       PetscInt m,n,M,N,mi,ni,Mi,Ni;
1313       Mat      B = vs->m[i][j];
1314       if (!B) continue;
1315       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
1316       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
1317       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
1318       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
1319       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
1320       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1321       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);
1322       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);
1323     }
1324   }
1325 #endif
1326 
1327   /* Set A->assembled if all non-null blocks are currently assembled */
1328   for (i=0; i<vs->nr; i++) {
1329     for (j=0; j<vs->nc; j++) {
1330       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1331     }
1332   }
1333   A->assembled = PETSC_TRUE;
1334   PetscFunctionReturn(0);
1335 }
1336 
1337 #undef __FUNCT__
1338 #define __FUNCT__ "MatCreateNest"
1339 /*@C
1340    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1341 
1342    Collective on Mat
1343 
1344    Input Parameter:
1345 +  comm - Communicator for the new Mat
1346 .  nr - number of nested row blocks
1347 .  is_row - index sets for each nested row block, or NULL to make contiguous
1348 .  nc - number of nested column blocks
1349 .  is_col - index sets for each nested column block, or NULL to make contiguous
1350 -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1351 
1352    Output Parameter:
1353 .  B - new matrix
1354 
1355    Level: advanced
1356 
1357 .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1358 @*/
1359 PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1360 {
1361   Mat            A;
1362   PetscErrorCode ierr;
1363 
1364   PetscFunctionBegin;
1365   *B   = 0;
1366   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1367   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1368   ierr = MatSetUp(A);CHKERRQ(ierr);
1369   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1370   *B   = A;
1371   PetscFunctionReturn(0);
1372 }
1373 
1374 #undef __FUNCT__
1375 #define __FUNCT__ "MatConvert_Nest_AIJ"
1376 PETSC_EXTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1377 {
1378   PetscErrorCode ierr;
1379   Mat_Nest       *nest = (Mat_Nest*)A->data;
1380   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz;
1381   Mat            C;
1382 
1383   PetscFunctionBegin;
1384   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1385   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1386   switch (reuse) {
1387   case MAT_INITIAL_MATRIX:
1388     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1389     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1390     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1391     *newmat = C;
1392     break;
1393   case MAT_REUSE_MATRIX:
1394     C = *newmat;
1395     break;
1396   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1397   }
1398 
1399   /* Preallocation */
1400   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1401   onnz = dnnz + m;
1402   for (k=0; k<m; k++) {
1403     dnnz[k] = 0;
1404     onnz[k] = 0;
1405   }
1406   for (j=0; j<nest->nc; ++j) {
1407     IS             bNis;
1408     PetscInt       bN;
1409     const PetscInt *bNindices;
1410     /* Using global column indices and ISAllGather() is not scalable. */
1411     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1412     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1413     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1414     for (i=0; i<nest->nr; ++i) {
1415       PetscSF        bmsf;
1416       PetscSFNode    *bmedges;
1417       Mat            B;
1418       PetscInt       bm, *bmdnnz, br;
1419       const PetscInt *bmindices;
1420       B = nest->m[i][j];
1421       if (!B) continue;
1422       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1423       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1424       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1425       ierr = PetscMalloc1(2*bm,&bmedges);CHKERRQ(ierr);
1426       ierr = PetscMalloc1(2*bm,&bmdnnz);CHKERRQ(ierr);
1427       for (k = 0; k < 2*bm; ++k) bmdnnz[k] = 0;
1428       /*
1429        Locate the owners for all of the locally-owned global row indices for this row block.
1430        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1431        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1432        */
1433       for (br = 0; br < bm; ++br) {
1434         PetscInt       row = bmindices[br], rowowner = 0, brncols, col, colowner = 0;
1435         const PetscInt *brcols;
1436         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1437         PetscInt       rowownerm; /* local row size on row's owning rank. */
1438 
1439         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1440         rowownerm = A->rmap->range[rowowner+1]-A->rmap->range[rowowner];
1441 
1442         bmedges[br].rank = rowowner; bmedges[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1443         bmedges[br].rank = rowowner; bmedges[br].index = rowrel+rowownerm; /* edge from bmonnz to onnz */
1444         /* Now actually compute the data -- bmdnnz and bmonnz by looking at the global columns in the br row of this block. */
1445         /* Note that this is not a pessimistic bound only because we assume the index sets embedding the blocks do not overlap. */
1446         ierr = MatGetRow(B,br,&brncols,&brcols,NULL);CHKERRQ(ierr);
1447         for (k=0; k<brncols; k++) {
1448           col  = bNindices[brcols[k]];
1449           ierr = PetscLayoutFindOwnerIndex(A->cmap,col,&colowner,NULL);CHKERRQ(ierr);
1450           if (colowner == rowowner) bmdnnz[br]++;
1451           else onnz[br]++;
1452         }
1453         ierr = MatRestoreRow(B,br,&brncols,&brcols,NULL);CHKERRQ(ierr);
1454       }
1455       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1456       /* bsf will have to take care of disposing of bedges. */
1457       ierr = PetscSFSetGraph(bmsf,m,2*bm,NULL,PETSC_COPY_VALUES,bmedges,PETSC_OWN_POINTER);CHKERRQ(ierr);
1458       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);CHKERRQ(ierr);
1459       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);CHKERRQ(ierr);
1460       ierr = PetscFree(bmdnnz);CHKERRQ(ierr);
1461       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1462     }
1463     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1464     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1465   }
1466   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1467   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1468   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1469 
1470   /* Fill by row */
1471   for (j=0; j<nest->nc; ++j) {
1472     /* Using global column indices and ISAllGather() is not scalable. */
1473     IS             bNis;
1474     PetscInt       bN;
1475     const PetscInt *bNindices;
1476     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1477     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1478     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1479     for (i=0; i<nest->nr; ++i) {
1480       Mat            B;
1481       PetscInt       bm, br;
1482       const PetscInt *bmindices;
1483       B = nest->m[i][j];
1484       if (!B) continue;
1485       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1486       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1487       for (br = 0; br < bm; ++br) {
1488         PetscInt          row = bmindices[br], brncols,  *cols;
1489         const PetscInt    *brcols;
1490         const PetscScalar *brcoldata;
1491         ierr = MatGetRow(B,br,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1492         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
1493         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1494         /*
1495          Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1496          Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1497          */
1498         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
1499         ierr = MatRestoreRow(B,br,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1500         ierr = PetscFree(cols);CHKERRQ(ierr);
1501       }
1502       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1503     }
1504     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1505     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1506   }
1507   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1508   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1509   PetscFunctionReturn(0);
1510 }
1511 
1512 /*MC
1513   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1514 
1515   Level: intermediate
1516 
1517   Notes:
1518   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1519   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1520   It is usually used with DMComposite and DMCreateMatrix()
1521 
1522 .seealso: MatCreate(), MatType, MatCreateNest()
1523 M*/
1524 #undef __FUNCT__
1525 #define __FUNCT__ "MatCreate_Nest"
1526 PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1527 {
1528   Mat_Nest       *s;
1529   PetscErrorCode ierr;
1530 
1531   PetscFunctionBegin;
1532   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1533   A->data = (void*)s;
1534 
1535   s->nr            = -1;
1536   s->nc            = -1;
1537   s->m             = NULL;
1538   s->splitassembly = PETSC_FALSE;
1539 
1540   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1541 
1542   A->ops->mult                  = MatMult_Nest;
1543   A->ops->multadd               = MatMultAdd_Nest;
1544   A->ops->multtranspose         = MatMultTranspose_Nest;
1545   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1546   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1547   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1548   A->ops->zeroentries           = MatZeroEntries_Nest;
1549   A->ops->duplicate             = MatDuplicate_Nest;
1550   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1551   A->ops->destroy               = MatDestroy_Nest;
1552   A->ops->view                  = MatView_Nest;
1553   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1554   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1555   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1556   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1557   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1558   A->ops->scale                 = MatScale_Nest;
1559   A->ops->shift                 = MatShift_Nest;
1560 
1561   A->spptr        = 0;
1562   A->assembled    = PETSC_FALSE;
1563 
1564   /* expose Nest api's */
1565   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1566   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1567   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1568   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);CHKERRQ(ierr);
1569   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);CHKERRQ(ierr);
1570   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1571   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1572   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1573 
1574   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1575   PetscFunctionReturn(0);
1576 }
1577