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