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