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