xref: /petsc/src/mat/impls/nest/matnest.c (revision 01a79839fc82a7dabb7a87cd2a8bb532c6bfa88d)
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       vs->m[i] = PETSC_NULL;
261     }
262     ierr = PetscFree(vs->m);CHKERRQ(ierr);
263     vs->m = PETSC_NULL;
264   }
265   ierr = PetscFree(vs);CHKERRQ(ierr);
266 
267   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "",0);CHKERRQ(ierr);
268   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
269   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
270   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
271   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
272   PetscFunctionReturn(0);
273 }
274 
275 #undef __FUNCT__
276 #define __FUNCT__ "MatAssemblyBegin_Nest"
277 static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
278 {
279   Mat_Nest       *vs = (Mat_Nest*)A->data;
280   PetscInt       i,j;
281   PetscErrorCode ierr;
282 
283   PetscFunctionBegin;
284   for (i=0; i<vs->nr; i++) {
285     for (j=0; j<vs->nc; j++) {
286       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
287     }
288   }
289   PetscFunctionReturn(0);
290 }
291 
292 #undef __FUNCT__
293 #define __FUNCT__ "MatAssemblyEnd_Nest"
294 static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
295 {
296   Mat_Nest       *vs = (Mat_Nest*)A->data;
297   PetscInt       i,j;
298   PetscErrorCode ierr;
299 
300   PetscFunctionBegin;
301   for (i=0; i<vs->nr; i++) {
302     for (j=0; j<vs->nc; j++) {
303       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
304     }
305   }
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
311 static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
312 {
313   PetscErrorCode ierr;
314   Mat_Nest       *vs = (Mat_Nest*)A->data;
315   PetscInt       j;
316   Mat            sub;
317 
318   PetscFunctionBegin;
319   sub = vs->m[row][row];        /* Prefer to find on the diagonal */
320   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
321   if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",row);
322   ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */
323   *B = sub;
324   PetscFunctionReturn(0);
325 }
326 
327 #undef __FUNCT__
328 #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
329 static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
330 {
331   PetscErrorCode ierr;
332   Mat_Nest       *vs = (Mat_Nest*)A->data;
333   PetscInt       i;
334   Mat            sub;
335 
336   PetscFunctionBegin;
337   sub = vs->m[col][col];        /* Prefer to find on the diagonal */
338   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
339   if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",col);
340   ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */
341   *B = sub;
342   PetscFunctionReturn(0);
343 }
344 
345 #undef __FUNCT__
346 #define __FUNCT__ "MatNestFindIS"
347 static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
348 {
349   PetscErrorCode ierr;
350   PetscInt       i;
351   PetscBool      flg;
352 
353   PetscFunctionBegin;
354   PetscValidPointer(list,3);
355   PetscValidHeaderSpecific(is,IS_CLASSID,4);
356   PetscValidIntPointer(found,5);
357   *found = -1;
358   for (i=0; i<n; i++) {
359     if (!list[i]) continue;
360     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
361     if (flg) {
362       *found = i;
363       PetscFunctionReturn(0);
364     }
365   }
366   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "MatNestFindSubMat"
372 static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
373 {
374   PetscErrorCode ierr;
375   Mat_Nest       *vs = (Mat_Nest*)A->data;
376   PetscInt       row,col;
377 
378   PetscFunctionBegin;
379   ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
380   ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
381   *B = vs->m[row][col];
382   PetscFunctionReturn(0);
383 }
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "MatGetSubMatrix_Nest"
387 static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
388 {
389   PetscErrorCode ierr;
390   Mat_Nest       *vs = (Mat_Nest*)A->data;
391   Mat            sub;
392 
393   PetscFunctionBegin;
394   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
395   switch (reuse) {
396   case MAT_INITIAL_MATRIX:
397     ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);
398     *B = sub;
399     break;
400   case MAT_REUSE_MATRIX:
401     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
402     break;
403   case MAT_IGNORE_MATRIX:       /* Nothing to do */
404     break;
405   }
406   PetscFunctionReturn(0);
407 }
408 
409 #undef __FUNCT__
410 #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
411 PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
412 {
413   PetscErrorCode ierr;
414   Mat_Nest       *vs = (Mat_Nest*)A->data;
415   Mat            sub;
416 
417   PetscFunctionBegin;
418   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
419   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
420   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
421   *B = sub;
422   PetscFunctionReturn(0);
423 }
424 
425 #undef __FUNCT__
426 #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
427 static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
428 {
429   PetscErrorCode ierr;
430   Mat_Nest       *vs = (Mat_Nest*)A->data;
431   Mat            sub;
432 
433   PetscFunctionBegin;
434   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
435   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
436   if (sub) {
437     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
438     ierr = MatDestroy(*B);CHKERRQ(ierr);
439   }
440   *B = PETSC_NULL;
441   PetscFunctionReturn(0);
442 }
443 
444 #undef __FUNCT__
445 #define __FUNCT__ "MatGetVecs_Nest"
446 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
447 {
448   Mat_Nest       *bA = (Mat_Nest*)A->data;
449   Vec            *L,*R;
450   MPI_Comm       comm;
451   PetscInt       i,j;
452   PetscErrorCode ierr;
453 
454   PetscFunctionBegin;
455   comm = ((PetscObject)A)->comm;
456   if (right) {
457     /* allocate R */
458     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
459     /* Create the right vectors */
460     for (j=0; j<bA->nc; j++) {
461       for (i=0; i<bA->nr; i++) {
462         if (bA->m[i][j]) {
463           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
464           break;
465         }
466       }
467       if (i==bA->nr) {
468         /* have an empty column */
469         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
470       }
471     }
472     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
473     /* hand back control to the nest vector */
474     for (j=0; j<bA->nc; j++) {
475       ierr = VecDestroy(R[j]);CHKERRQ(ierr);
476     }
477     ierr = PetscFree(R);CHKERRQ(ierr);
478   }
479 
480   if (left) {
481     /* allocate L */
482     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
483     /* Create the left vectors */
484     for (i=0; i<bA->nr; i++) {
485       for (j=0; j<bA->nc; j++) {
486         if (bA->m[i][j]) {
487           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
488           break;
489         }
490       }
491       if (j==bA->nc) {
492         /* have an empty row */
493         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
494       }
495     }
496 
497     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
498     for (i=0; i<bA->nr; i++) {
499       ierr = VecDestroy(L[i]);CHKERRQ(ierr);
500     }
501 
502     ierr = PetscFree(L);CHKERRQ(ierr);
503   }
504   PetscFunctionReturn(0);
505 }
506 
507 #undef __FUNCT__
508 #define __FUNCT__ "MatView_Nest"
509 static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
510 {
511   Mat_Nest       *bA = (Mat_Nest*)A->data;
512   PetscBool      isascii;
513   PetscInt       i,j;
514   PetscErrorCode ierr;
515 
516   PetscFunctionBegin;
517   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
518   if (isascii) {
519 
520     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
521     PetscViewerASCIIPushTab( viewer );    /* push0 */
522     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
523 
524     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
525     for (i=0; i<bA->nr; i++) {
526       for (j=0; j<bA->nc; j++) {
527         const MatType type;
528         const char *name;
529         PetscInt NR,NC;
530         PetscBool isNest = PETSC_FALSE;
531 
532         if (!bA->m[i][j]) {
533           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
534           continue;
535         }
536         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
537         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
538         name = ((PetscObject)bA->m[i][j])->prefix;
539         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
540 
541         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
542 
543         if (isNest) {
544           PetscViewerASCIIPushTab( viewer );  /* push1 */
545           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
546           PetscViewerASCIIPopTab(viewer);    /* pop1 */
547         }
548       }
549     }
550     PetscViewerASCIIPopTab(viewer);    /* pop0 */
551   }
552   PetscFunctionReturn(0);
553 }
554 
555 #undef __FUNCT__
556 #define __FUNCT__ "MatZeroEntries_Nest"
557 static PetscErrorCode MatZeroEntries_Nest(Mat A)
558 {
559   Mat_Nest       *bA = (Mat_Nest*)A->data;
560   PetscInt       i,j;
561   PetscErrorCode ierr;
562 
563   PetscFunctionBegin;
564   for (i=0; i<bA->nr; i++) {
565     for (j=0; j<bA->nc; j++) {
566       if (!bA->m[i][j]) continue;
567       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
568     }
569   }
570   PetscFunctionReturn(0);
571 }
572 
573 #undef __FUNCT__
574 #define __FUNCT__ "MatDuplicate_Nest"
575 static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
576 {
577   Mat_Nest       *bA = (Mat_Nest*)A->data;
578   Mat            *b;
579   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
580   PetscErrorCode ierr;
581 
582   PetscFunctionBegin;
583   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
584   for (i=0; i<nr; i++) {
585     for (j=0; j<nc; j++) {
586       if (bA->m[i][j]) {
587         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
588       } else {
589         b[i*nc+j] = PETSC_NULL;
590       }
591     }
592   }
593   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
594   /* Give the new MatNest exclusive ownership */
595   for (i=0; i<nr*nc; i++) {
596     if (b[i]) {ierr = MatDestroy(b[i]);CHKERRQ(ierr);}
597   }
598   ierr = PetscFree(b);CHKERRQ(ierr);
599 
600   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
601   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
602   PetscFunctionReturn(0);
603 }
604 
605 /* nest api */
606 EXTERN_C_BEGIN
607 #undef __FUNCT__
608 #define __FUNCT__ "MatNestGetSubMat_Nest"
609 PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
610 {
611   Mat_Nest *bA = (Mat_Nest*)A->data;
612   PetscFunctionBegin;
613   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
614   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
615   *mat = bA->m[idxm][jdxm];
616   PetscFunctionReturn(0);
617 }
618 EXTERN_C_END
619 
620 #undef __FUNCT__
621 #define __FUNCT__ "MatNestGetSubMat"
622 /*@C
623  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
624 
625  Not collective
626 
627  Input Parameters:
628  .  A  - nest matrix
629  .  idxm - index of the matrix within the nest matrix
630  .  jdxm - index of the matrix within the nest matrix
631 
632  Output Parameter:
633  .  sub - matrix at index idxm,jdxm within the nest matrix
634 
635  Notes:
636 
637  Level: developer
638 
639  .seealso: MatNestGetSize(), MatNestGetSubMats()
640 @*/
641 PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
642 {
643   PetscErrorCode ierr;
644 
645   PetscFunctionBegin;
646   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
647   PetscFunctionReturn(0);
648 }
649 
650 EXTERN_C_BEGIN
651 #undef __FUNCT__
652 #define __FUNCT__ "MatNestGetSubMats_Nest"
653 PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
654 {
655   Mat_Nest *bA = (Mat_Nest*)A->data;
656   PetscFunctionBegin;
657   if (M)   { *M   = bA->nr; }
658   if (N)   { *N   = bA->nc; }
659   if (mat) { *mat = bA->m;  }
660   PetscFunctionReturn(0);
661 }
662 EXTERN_C_END
663 
664 #undef __FUNCT__
665 #define __FUNCT__ "MatNestGetSubMats"
666 /*@C
667  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
668 
669  Not collective
670 
671  Input Parameters:
672  .  A  - nest matri
673 
674  Output Parameter:
675  .  M - number of rows in the nest matrix
676  .  N - number of cols in the nest matrix
677  .  mat - 2d array of matrices
678 
679  Notes:
680 
681  The user should not free the array mat.
682 
683  Level: developer
684 
685  .seealso: MatNestGetSize(), MatNestGetSubMat()
686 @*/
687 PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
688 {
689   PetscErrorCode ierr;
690 
691   PetscFunctionBegin;
692   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
693   PetscFunctionReturn(0);
694 }
695 
696 EXTERN_C_BEGIN
697 #undef __FUNCT__
698 #define __FUNCT__ "MatNestGetSize_Nest"
699 PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
700 {
701   Mat_Nest  *bA = (Mat_Nest*)A->data;
702 
703   PetscFunctionBegin;
704   if (M) { *M  = bA->nr; }
705   if (N) { *N  = bA->nc; }
706   PetscFunctionReturn(0);
707 }
708 EXTERN_C_END
709 
710 #undef __FUNCT__
711 #define __FUNCT__ "MatNestGetSize"
712 /*@C
713  MatNestGetSize - Returns the size of the nest matrix.
714 
715  Not collective
716 
717  Input Parameters:
718  .  A  - nest matrix
719 
720  Output Parameter:
721  .  M - number of rows in the nested mat
722  .  N - number of cols in the nested mat
723 
724  Notes:
725 
726  Level: developer
727 
728  .seealso: MatNestGetSubMat(), MatNestGetSubMats()
729 @*/
730 PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
731 {
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
736   PetscFunctionReturn(0);
737 }
738 
739 EXTERN_C_BEGIN
740 #undef __FUNCT__
741 #define __FUNCT__ "MatNestSetVecType_Nest"
742 PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
743 {
744   PetscErrorCode ierr;
745   PetscBool      flg;
746 
747   PetscFunctionBegin;
748   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
749   /* In reality, this only distinguishes VECNEST and "other" */
750   A->ops->getvecs = flg ? MatGetVecs_Nest : 0;
751   PetscFunctionReturn(0);
752 }
753 EXTERN_C_END
754 
755 #undef __FUNCT__
756 #define __FUNCT__ "MatNestSetVecType"
757 /*@C
758  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
759 
760  Not collective
761 
762  Input Parameters:
763 +  A  - nest matrix
764 -  vtype - type to use for creating vectors
765 
766  Notes:
767 
768  Level: developer
769 
770  .seealso: MatGetVecs()
771 @*/
772 PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
773 {
774   PetscErrorCode ierr;
775 
776   PetscFunctionBegin;
777   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
778   PetscFunctionReturn(0);
779 }
780 
781 EXTERN_C_BEGIN
782 #undef __FUNCT__
783 #define __FUNCT__ "MatNestSetSubMats_Nest"
784 PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
785 {
786   Mat_Nest       *s = (Mat_Nest*)A->data;
787   PetscInt       i,j,m,n,M,N;
788   PetscErrorCode ierr;
789 
790   PetscFunctionBegin;
791   s->nr = nr;
792   s->nc = nc;
793 
794   /* Create space for submatrices */
795   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
796   for (i=0; i<nr; i++) {
797     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
798   }
799   for (i=0; i<nr; i++) {
800     for (j=0; j<nc; j++) {
801       s->m[i][j] = a[i*nc+j];
802       if (a[i*nc+j]) {
803         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
804       }
805     }
806   }
807 
808   ierr = PetscMalloc(sizeof(IS)*nr,&s->isglobal.row);CHKERRQ(ierr);
809   ierr = PetscMalloc(sizeof(IS)*nc,&s->isglobal.col);CHKERRQ(ierr);
810 
811   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
812   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
813   for (i=0; i<nr; i++) s->row_len[i]=-1;
814   for (j=0; j<nc; j++) s->col_len[j]=-1;
815 
816   m = n = 0;
817   M = N = 0;
818   ierr = MatNestGetSize_Private(A,PETSC_TRUE,&M,&N);CHKERRQ(ierr);
819   ierr = MatNestGetSize_Private(A,PETSC_FALSE,&m,&n);CHKERRQ(ierr);
820 
821   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
822   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
823   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
824   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
825   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
826   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
827 
828   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
829   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
830 
831   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
832   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
833   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
834   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
835   PetscFunctionReturn(0);
836 }
837 EXTERN_C_END
838 
839 #undef __FUNCT__
840 #define __FUNCT__ "MatNestSetSubMats"
841 /*@
842    MatNestSetSubMats - Sets the nested submatrices
843 
844    Collective on Mat
845 
846    Input Parameter:
847 +  N - nested matrix
848 .  nr - number of nested row blocks
849 .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
850 .  nc - number of nested column blocks
851 .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
852 -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
853 
854    Level: advanced
855 
856 .seealso: MatCreateNest(), MATNEST
857 @*/
858 PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
859 {
860   PetscErrorCode ierr;
861   PetscInt i;
862 
863   PetscFunctionBegin;
864   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
865   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
866   if (nr && is_row) {
867     PetscValidPointer(is_row,3);
868     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
869   }
870   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
871   if (nc && is_row) {
872     PetscValidPointer(is_col,5);
873     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
874   }
875   if (nr*nc) PetscValidPointer(a,6);
876   ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr);
877   PetscFunctionReturn(0);
878 }
879 
880 /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
881 /*
882   nprocessors = NP
883   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
884        proc 0: => (g_0,h_0,)
885        proc 1: => (g_1,h_1,)
886        ...
887        proc nprocs-1: => (g_NP-1,h_NP-1,)
888 
889             proc 0:                      proc 1:                    proc nprocs-1:
890     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
891 
892             proc 0:
893     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
894             proc 1:
895     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
896 
897             proc NP-1:
898     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
899 */
900 #undef __FUNCT__
901 #define __FUNCT__ "MatSetUp_NestIS_Private"
902 static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
903 {
904   Mat_Nest       *vs = (Mat_Nest*)A->data;
905   PetscInt       i,j,offset,n,bs;
906   PetscErrorCode ierr;
907   Mat            sub;
908 
909   PetscFunctionBegin;
910   if (is_row) { /* valid IS is passed in */
911     /* refs on is[] are incremeneted */
912     for (i=0; i<vs->nr; i++) {
913       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
914       vs->isglobal.row[i] = is_row[i];
915     }
916   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
917     offset = A->rmap->rstart;
918     for (i=0; i<vs->nr; i++) {
919       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
920       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
921       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
922       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
923       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
924       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
925       offset += n;
926     }
927   }
928 
929   if (is_col) { /* valid IS is passed in */
930     /* refs on is[] are incremeneted */
931     for (j=0; j<vs->nc; j++) {
932       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
933       vs->isglobal.col[j] = is_col[j];
934     }
935   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
936     offset = A->cmap->rstart;
937     for (j=0; j<vs->nc; j++) {
938       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
939       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
940       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
941       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
942       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
943       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
944       offset += n;
945     }
946   }
947 
948   /* Set up the local ISs */
949   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
950   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
951   for (i=0,offset=0; i<vs->nr; i++) {
952     IS                     isloc;
953     ISLocalToGlobalMapping rmap;
954     PetscInt               nlocal,bs;
955     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
956     ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);
957     if (rmap) {
958       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
959       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
960       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
961       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
962     } else {
963       nlocal = 0;
964       isloc  = PETSC_NULL;
965     }
966     vs->islocal.row[i] = isloc;
967     offset += nlocal;
968   }
969   for (i=0,offset=0; i<vs->nr; i++) {
970     IS                     isloc;
971     ISLocalToGlobalMapping cmap;
972     PetscInt               nlocal,bs;
973     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
974     ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);
975     if (cmap) {
976       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
977       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
978       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
979       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
980     } else {
981       nlocal = 0;
982       isloc  = PETSC_NULL;
983     }
984     vs->islocal.col[i] = isloc;
985     offset += nlocal;
986   }
987   PetscFunctionReturn(0);
988 }
989 
990 #undef __FUNCT__
991 #define __FUNCT__ "MatCreateNest"
992 /*@
993    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
994 
995    Collective on Mat
996 
997    Input Parameter:
998 +  comm - Communicator for the new Mat
999 .  nr - number of nested row blocks
1000 .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1001 .  nc - number of nested column blocks
1002 .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1003 -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1004 
1005    Output Parameter:
1006 .  B - new matrix
1007 
1008    Level: advanced
1009 
1010 .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST
1011 @*/
1012 PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1013 {
1014   Mat            A;
1015   PetscErrorCode ierr;
1016 
1017   PetscFunctionBegin;
1018   *B = 0;
1019   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1020   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1021   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1022   *B = A;
1023   PetscFunctionReturn(0);
1024 }
1025 
1026 /*MC
1027   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1028 
1029   Level: intermediate
1030 
1031   Notes:
1032   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1033   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1034   It is usually used with DMComposite and DMGetMatrix()
1035 
1036 .seealso: MatCreate(), MatType, MatCreateNest()
1037 M*/
1038 EXTERN_C_BEGIN
1039 #undef __FUNCT__
1040 #define __FUNCT__ "MatCreate_Nest"
1041 PetscErrorCode MatCreate_Nest(Mat A)
1042 {
1043   Mat_Nest       *s;
1044   PetscErrorCode ierr;
1045 
1046   PetscFunctionBegin;
1047   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1048   A->data         = (void*)s;
1049   s->nr = s->nc   = -1;
1050   s->m            = PETSC_NULL;
1051 
1052   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1053   A->ops->mult                  = MatMult_Nest;
1054   A->ops->multadd               = 0;
1055   A->ops->multtranspose         = MatMultTranspose_Nest;
1056   A->ops->multtransposeadd      = 0;
1057   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1058   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1059   A->ops->zeroentries           = MatZeroEntries_Nest;
1060   A->ops->duplicate             = MatDuplicate_Nest;
1061   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1062   A->ops->destroy               = MatDestroy_Nest;
1063   A->ops->view                  = MatView_Nest;
1064   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1065   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1066   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1067 
1068   A->spptr        = 0;
1069   A->same_nonzero = PETSC_FALSE;
1070   A->assembled    = PETSC_FALSE;
1071 
1072   /* expose Nest api's */
1073   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1074   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1075   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1076   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1077   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1078 
1079   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1080   PetscFunctionReturn(0);
1081 }
1082 EXTERN_C_END
1083