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