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