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