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