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