xref: /petsc/src/mat/impls/nest/matnest.c (revision e00f73327efb23b4e84a88ac93064d8a8f58c4e4)
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,"MatNestSetSubMat_C",   "",0);CHKERRQ(ierr);
254   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
255   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
256   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
257   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "MatAssemblyBegin_Nest"
263 static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
264 {
265   Mat_Nest       *vs = (Mat_Nest*)A->data;
266   PetscInt       i,j;
267   PetscErrorCode ierr;
268 
269   PetscFunctionBegin;
270   for (i=0; i<vs->nr; i++) {
271     for (j=0; j<vs->nc; j++) {
272       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
273     }
274   }
275   PetscFunctionReturn(0);
276 }
277 
278 #undef __FUNCT__
279 #define __FUNCT__ "MatAssemblyEnd_Nest"
280 static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
281 {
282   Mat_Nest       *vs = (Mat_Nest*)A->data;
283   PetscInt       i,j;
284   PetscErrorCode ierr;
285 
286   PetscFunctionBegin;
287   for (i=0; i<vs->nr; i++) {
288     for (j=0; j<vs->nc; j++) {
289       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
290     }
291   }
292   PetscFunctionReturn(0);
293 }
294 
295 #undef __FUNCT__
296 #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
297 static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
298 {
299   PetscErrorCode ierr;
300   Mat_Nest       *vs = (Mat_Nest*)A->data;
301   PetscInt       j;
302   Mat            sub;
303 
304   PetscFunctionBegin;
305   sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */
306   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
307   if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */
308   *B = sub;
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
314 static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
315 {
316   PetscErrorCode ierr;
317   Mat_Nest       *vs = (Mat_Nest*)A->data;
318   PetscInt       i;
319   Mat            sub;
320 
321   PetscFunctionBegin;
322   sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */
323   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
324   if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */
325   *B = sub;
326   PetscFunctionReturn(0);
327 }
328 
329 #undef __FUNCT__
330 #define __FUNCT__ "MatNestFindIS"
331 static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
332 {
333   PetscErrorCode ierr;
334   PetscInt       i;
335   PetscBool      flg;
336 
337   PetscFunctionBegin;
338   PetscValidPointer(list,3);
339   PetscValidHeaderSpecific(is,IS_CLASSID,4);
340   PetscValidIntPointer(found,5);
341   *found = -1;
342   for (i=0; i<n; i++) {
343     if (!list[i]) continue;
344     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
345     if (flg) {
346       *found = i;
347       PetscFunctionReturn(0);
348     }
349   }
350   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
351   PetscFunctionReturn(0);
352 }
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "MatNestGetRow"
356 /* Get a block row as a new MatNest */
357 static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
358 {
359   Mat_Nest       *vs = (Mat_Nest*)A->data;
360   Mat            C;
361   char           keyname[256];
362   PetscErrorCode ierr;
363 
364   PetscFunctionBegin;
365   *B = PETSC_NULL;
366   ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr);
367   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
368   if (*B) PetscFunctionReturn(0);
369 
370   ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
371   (*B)->assembled = A->assembled;
372   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
373   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
374   PetscFunctionReturn(0);
375 }
376 
377 #undef __FUNCT__
378 #define __FUNCT__ "MatNestFindSubMat"
379 static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
380 {
381   Mat_Nest       *vs = (Mat_Nest*)A->data;
382   PetscErrorCode ierr;
383   PetscInt       row,col,i;
384   IS             *iscopy;
385   Mat            *matcopy;
386   PetscBool      same,isFullCol;
387 
388   PetscFunctionBegin;
389   /* Check if full column space. This is a hack */
390   isFullCol = PETSC_FALSE;
391   ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
392   if (same) {
393     PetscInt n,first,step;
394     ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
395     ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
396     if (A->cmap->n == n && A->cmap->rstart == first && step == 1) isFullCol = PETSC_TRUE;
397   }
398 
399   if (isFullCol) {
400     PetscInt row;
401     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
402     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
403   } else {
404     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
405     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
406     *B = vs->m[row][col];
407   }
408   PetscFunctionReturn(0);
409 }
410 
411 #undef __FUNCT__
412 #define __FUNCT__ "MatGetSubMatrix_Nest"
413 static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
414 {
415   PetscErrorCode ierr;
416   Mat_Nest       *vs = (Mat_Nest*)A->data;
417   Mat            sub;
418 
419   PetscFunctionBegin;
420   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
421   switch (reuse) {
422   case MAT_INITIAL_MATRIX:
423     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
424     *B = sub;
425     break;
426   case MAT_REUSE_MATRIX:
427     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
428     break;
429   case MAT_IGNORE_MATRIX:       /* Nothing to do */
430     break;
431   }
432   PetscFunctionReturn(0);
433 }
434 
435 #undef __FUNCT__
436 #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
437 PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
438 {
439   PetscErrorCode ierr;
440   Mat_Nest       *vs = (Mat_Nest*)A->data;
441   Mat            sub;
442 
443   PetscFunctionBegin;
444   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
445   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
446   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
447   *B = sub;
448   PetscFunctionReturn(0);
449 }
450 
451 #undef __FUNCT__
452 #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
453 static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
454 {
455   PetscErrorCode ierr;
456   Mat_Nest       *vs = (Mat_Nest*)A->data;
457   Mat            sub;
458 
459   PetscFunctionBegin;
460   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
461   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
462   if (sub) {
463     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
464     ierr = MatDestroy(B);CHKERRQ(ierr);
465   }
466   PetscFunctionReturn(0);
467 }
468 
469 #undef __FUNCT__
470 #define __FUNCT__ "MatGetDiagonal_Nest"
471 static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
472 {
473   Mat_Nest       *bA = (Mat_Nest*)A->data;
474   Vec            *bdiag;
475   PetscInt       i;
476   PetscErrorCode ierr;
477 
478   PetscFunctionBegin;
479   /*  ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr); */
480   /*  ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr); */
481   ierr = VecNestGetSubVecs(v,PETSC_NULL,&bdiag);CHKERRQ(ierr);
482   for (i=0; i<bA->nr; i++) {
483     if (bA->m[i][i]) {
484       ierr = MatGetDiagonal(bA->m[i][i],bdiag[i]);CHKERRQ(ierr);
485     } else {
486       ierr = VecSet(bdiag[i],1.0);CHKERRQ(ierr);
487     }
488   }
489   PetscFunctionReturn(0);
490 }
491 
492 #undef __FUNCT__
493 #define __FUNCT__ "MatGetDiagonal_Nest2"
494 static PetscErrorCode MatGetDiagonal_Nest2(Mat A,Vec v)
495 {
496   Mat_Nest       *bA = (Mat_Nest*)A->data;
497   Vec            diag,*bdiag;
498   VecScatter     *vscat;
499   PetscInt       i;
500   PetscErrorCode ierr;
501 
502   PetscFunctionBegin;
503   ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr);
504   ierr = MatGetDiagonal_Nest(A,diag);CHKERRQ(ierr);
505 
506   /* scatter diag into v */
507   ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr);
508   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscat );CHKERRQ(ierr);
509   for (i=0; i<bA->nr; i++) {
510     ierr = VecScatterCreate(v,bA->isglobal.row[i], bdiag[i],PETSC_NULL,&vscat[i]);CHKERRQ(ierr);
511     ierr = VecScatterBegin(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
512   }
513   for (i=0; i<bA->nr; i++) {
514     ierr = VecScatterEnd(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
515   }
516   for (i=0; i<bA->nr; i++) {
517     ierr = VecScatterDestroy(&vscat[i]);CHKERRQ(ierr);
518   }
519   ierr = PetscFree(vscat);CHKERRQ(ierr);
520   ierr = VecDestroy(&diag);CHKERRQ(ierr);
521   PetscFunctionReturn(0);
522 }
523 
524 #undef __FUNCT__
525 #define __FUNCT__ "MatDiagonalScale_Nest"
526 static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
527 {
528   Mat_Nest       *bA = (Mat_Nest*)A->data;
529   Vec            *bl,*br;
530   PetscInt       i,j;
531   PetscErrorCode ierr;
532 
533   PetscFunctionBegin;
534   ierr = VecNestGetSubVecs(l,PETSC_NULL,&bl);CHKERRQ(ierr);
535   ierr = VecNestGetSubVecs(r,PETSC_NULL,&br);CHKERRQ(ierr);
536   for (i=0; i<bA->nr; i++) {
537     for (j=0; j<bA->nc; j++) {
538       if (bA->m[i][j]) {
539         ierr = MatDiagonalScale(bA->m[i][j],bl[i],br[j]);CHKERRQ(ierr);
540       }
541     }
542   }
543   PetscFunctionReturn(0);
544 }
545 
546 #undef __FUNCT__
547 #define __FUNCT__ "MatDiagonalScale_Nest2"
548 static PetscErrorCode MatDiagonalScale_Nest2(Mat A,Vec l,Vec r)
549 {
550   Mat_Nest       *bA = (Mat_Nest*)A->data;
551   Vec            bl,br,*ble,*bre;
552   VecScatter     *vscatl,*vscatr;
553   PetscInt       i;
554   PetscErrorCode ierr;
555 
556   PetscFunctionBegin;
557   /* scatter l,r into bl,br */
558   ierr = MatGetVecs_Nest(A,&bl,&br);CHKERRQ(ierr);
559   ierr = VecNestGetSubVecs(bl,PETSC_NULL,&ble);CHKERRQ(ierr);
560   ierr = VecNestGetSubVecs(br,PETSC_NULL,&bre);CHKERRQ(ierr);
561 
562   /* row */
563   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscatl );CHKERRQ(ierr);
564   for (i=0; i<bA->nr; i++) {
565     ierr = VecScatterCreate(l,bA->isglobal.row[i], ble[i],PETSC_NULL,&vscatl[i]);CHKERRQ(ierr);
566     ierr = VecScatterBegin(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
567   }
568   for (i=0; i<bA->nr; i++) {
569     ierr = VecScatterEnd(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
570   }
571   /* col */
572   ierr = PetscMalloc( sizeof(VecScatter) * bA->nc, &vscatr );CHKERRQ(ierr);
573   for (i=0; i<bA->nc; i++) {
574     ierr = VecScatterCreate(r,bA->isglobal.col[i], bre[i],PETSC_NULL,&vscatr[i]);CHKERRQ(ierr);
575     ierr = VecScatterBegin(vscatr[i],l,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
576   }
577   for (i=0; i<bA->nc; i++) {
578     ierr = VecScatterEnd(vscatr[i],r,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
579   }
580 
581   ierr = MatDiagonalScale_Nest(A,bl,br);CHKERRQ(ierr);
582 
583   for (i=0; i<bA->nr; i++) {
584     ierr = VecScatterDestroy(&vscatl[i]);CHKERRQ(ierr);
585   }
586   for (i=0; i<bA->nc; i++) {
587     ierr = VecScatterDestroy(&vscatr[i]);CHKERRQ(ierr);
588   }
589   ierr = PetscFree(vscatl);CHKERRQ(ierr);
590   ierr = PetscFree(vscatr);CHKERRQ(ierr);
591   ierr = VecDestroy(&bl);CHKERRQ(ierr);
592   ierr = VecDestroy(&br);CHKERRQ(ierr);
593   PetscFunctionReturn(0);
594 }
595 
596 #undef __FUNCT__
597 #define __FUNCT__ "MatGetVecs_Nest"
598 static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
599 {
600   Mat_Nest       *bA = (Mat_Nest*)A->data;
601   Vec            *L,*R;
602   MPI_Comm       comm;
603   PetscInt       i,j;
604   PetscErrorCode ierr;
605 
606   PetscFunctionBegin;
607   comm = ((PetscObject)A)->comm;
608   if (right) {
609     /* allocate R */
610     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
611     /* Create the right vectors */
612     for (j=0; j<bA->nc; j++) {
613       for (i=0; i<bA->nr; i++) {
614         if (bA->m[i][j]) {
615           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
616           break;
617         }
618       }
619       if (i==bA->nr) {
620         /* have an empty column */
621         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
622       }
623     }
624     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
625     /* hand back control to the nest vector */
626     for (j=0; j<bA->nc; j++) {
627       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
628     }
629     ierr = PetscFree(R);CHKERRQ(ierr);
630   }
631 
632   if (left) {
633     /* allocate L */
634     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
635     /* Create the left vectors */
636     for (i=0; i<bA->nr; i++) {
637       for (j=0; j<bA->nc; j++) {
638         if (bA->m[i][j]) {
639           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
640           break;
641         }
642       }
643       if (j==bA->nc) {
644         /* have an empty row */
645         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
646       }
647     }
648 
649     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
650     for (i=0; i<bA->nr; i++) {
651       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
652     }
653 
654     ierr = PetscFree(L);CHKERRQ(ierr);
655   }
656   PetscFunctionReturn(0);
657 }
658 
659 #undef __FUNCT__
660 #define __FUNCT__ "MatView_Nest"
661 static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
662 {
663   Mat_Nest       *bA = (Mat_Nest*)A->data;
664   PetscBool      isascii;
665   PetscInt       i,j;
666   PetscErrorCode ierr;
667 
668   PetscFunctionBegin;
669   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
670   if (isascii) {
671 
672     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
673     PetscViewerASCIIPushTab( viewer );    /* push0 */
674     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
675 
676     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
677     for (i=0; i<bA->nr; i++) {
678       for (j=0; j<bA->nc; j++) {
679         const MatType type;
680         const char *name;
681         PetscInt NR,NC;
682         PetscBool isNest = PETSC_FALSE;
683 
684         if (!bA->m[i][j]) {
685           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
686           continue;
687         }
688         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
689         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
690         name = ((PetscObject)bA->m[i][j])->prefix;
691         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
692 
693         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
694 
695         if (isNest) {
696           PetscViewerASCIIPushTab( viewer );  /* push1 */
697           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
698           PetscViewerASCIIPopTab(viewer);    /* pop1 */
699         }
700       }
701     }
702     PetscViewerASCIIPopTab(viewer);    /* pop0 */
703   }
704   PetscFunctionReturn(0);
705 }
706 
707 #undef __FUNCT__
708 #define __FUNCT__ "MatZeroEntries_Nest"
709 static PetscErrorCode MatZeroEntries_Nest(Mat A)
710 {
711   Mat_Nest       *bA = (Mat_Nest*)A->data;
712   PetscInt       i,j;
713   PetscErrorCode ierr;
714 
715   PetscFunctionBegin;
716   for (i=0; i<bA->nr; i++) {
717     for (j=0; j<bA->nc; j++) {
718       if (!bA->m[i][j]) continue;
719       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
720     }
721   }
722   PetscFunctionReturn(0);
723 }
724 
725 #undef __FUNCT__
726 #define __FUNCT__ "MatDuplicate_Nest"
727 static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
728 {
729   Mat_Nest       *bA = (Mat_Nest*)A->data;
730   Mat            *b;
731   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
736   for (i=0; i<nr; i++) {
737     for (j=0; j<nc; j++) {
738       if (bA->m[i][j]) {
739         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
740       } else {
741         b[i*nc+j] = PETSC_NULL;
742       }
743     }
744   }
745   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
746   /* Give the new MatNest exclusive ownership */
747   for (i=0; i<nr*nc; i++) {
748     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
749   }
750   ierr = PetscFree(b);CHKERRQ(ierr);
751 
752   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
753   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
754   PetscFunctionReturn(0);
755 }
756 
757 /* nest api */
758 EXTERN_C_BEGIN
759 #undef __FUNCT__
760 #define __FUNCT__ "MatNestGetSubMat_Nest"
761 PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
762 {
763   Mat_Nest *bA = (Mat_Nest*)A->data;
764   PetscFunctionBegin;
765   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
766   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
767   *mat = bA->m[idxm][jdxm];
768   PetscFunctionReturn(0);
769 }
770 EXTERN_C_END
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "MatNestGetSubMat"
774 /*@C
775  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
776 
777  Not collective
778 
779  Input Parameters:
780 +   A  - nest matrix
781 .   idxm - index of the matrix within the nest matrix
782 -   jdxm - index of the matrix within the nest matrix
783 
784  Output Parameter:
785 .   sub - matrix at index idxm,jdxm within the nest matrix
786 
787  Level: developer
788 
789 .seealso: MatNestGetSize(), MatNestGetSubMats()
790 @*/
791 PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
792 {
793   PetscErrorCode ierr;
794 
795   PetscFunctionBegin;
796   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
797   PetscFunctionReturn(0);
798 }
799 
800 EXTERN_C_BEGIN
801 #undef __FUNCT__
802 #define __FUNCT__ "MatNestSetSubMat_Nest"
803 PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
804 {
805   Mat_Nest *bA = (Mat_Nest*)A->data;
806   PetscInt m,n,M,N,mi,ni,Mi,Ni;
807   PetscErrorCode ierr;
808 
809   PetscFunctionBegin;
810   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
811   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
812   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
813   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
814   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
815   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
816   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
817   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
818   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);
819   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);
820   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
821   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
822   bA->m[idxm][jdxm] = mat;
823   PetscFunctionReturn(0);
824 }
825 EXTERN_C_END
826 
827 #undef __FUNCT__
828 #define __FUNCT__ "MatNestSetSubMat"
829 /*@C
830  MatNestSetSubMat - Set a single submatrix in the nest matrix.
831 
832  Logically collective on the submatrix communicator
833 
834  Input Parameters:
835 +   A  - nest matrix
836 .   idxm - index of the matrix within the nest matrix
837 .   jdxm - index of the matrix within the nest matrix
838 -   sub - matrix at index idxm,jdxm within the nest matrix
839 
840  Notes:
841  The new submatrix must have the same size and communicator as that block of the nest.
842 
843  This increments the reference count of the submatrix.
844 
845  Level: developer
846 
847 .seealso: MatNestSetSubMats(), MatNestGetSubMat()
848 @*/
849 PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
850 {
851   PetscErrorCode ierr;
852 
853   PetscFunctionBegin;
854   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
855   PetscFunctionReturn(0);
856 }
857 
858 EXTERN_C_BEGIN
859 #undef __FUNCT__
860 #define __FUNCT__ "MatNestGetSubMats_Nest"
861 PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
862 {
863   Mat_Nest *bA = (Mat_Nest*)A->data;
864   PetscFunctionBegin;
865   if (M)   { *M   = bA->nr; }
866   if (N)   { *N   = bA->nc; }
867   if (mat) { *mat = bA->m;  }
868   PetscFunctionReturn(0);
869 }
870 EXTERN_C_END
871 
872 #undef __FUNCT__
873 #define __FUNCT__ "MatNestGetSubMats"
874 /*@C
875  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
876 
877  Not collective
878 
879  Input Parameters:
880 .   A  - nest matrix
881 
882  Output Parameter:
883 +   M - number of rows in the nest matrix
884 .   N - number of cols in the nest matrix
885 -   mat - 2d array of matrices
886 
887  Notes:
888 
889  The user should not free the array mat.
890 
891  Level: developer
892 
893 .seealso: MatNestGetSize(), MatNestGetSubMat()
894 @*/
895 PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
896 {
897   PetscErrorCode ierr;
898 
899   PetscFunctionBegin;
900   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
901   PetscFunctionReturn(0);
902 }
903 
904 EXTERN_C_BEGIN
905 #undef __FUNCT__
906 #define __FUNCT__ "MatNestGetSize_Nest"
907 PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
908 {
909   Mat_Nest  *bA = (Mat_Nest*)A->data;
910 
911   PetscFunctionBegin;
912   if (M) { *M  = bA->nr; }
913   if (N) { *N  = bA->nc; }
914   PetscFunctionReturn(0);
915 }
916 EXTERN_C_END
917 
918 #undef __FUNCT__
919 #define __FUNCT__ "MatNestGetSize"
920 /*@C
921  MatNestGetSize - Returns the size of the nest matrix.
922 
923  Not collective
924 
925  Input Parameters:
926 .   A  - nest matrix
927 
928  Output Parameter:
929 +   M - number of rows in the nested mat
930 -   N - number of cols in the nested mat
931 
932  Notes:
933 
934  Level: developer
935 
936 .seealso: MatNestGetSubMat(), MatNestGetSubMats()
937 @*/
938 PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
939 {
940   PetscErrorCode ierr;
941 
942   PetscFunctionBegin;
943   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
944   PetscFunctionReturn(0);
945 }
946 
947 EXTERN_C_BEGIN
948 #undef __FUNCT__
949 #define __FUNCT__ "MatNestSetVecType_Nest"
950 PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
951 {
952   PetscErrorCode ierr;
953   PetscBool      flg;
954 
955   PetscFunctionBegin;
956   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
957   /* In reality, this only distinguishes VECNEST and "other" */
958   A->ops->getvecs       = flg ? MatGetVecs_Nest       : 0;
959   A->ops->getdiagonal   = flg ? MatGetDiagonal_Nest   : 0;
960   A->ops->diagonalscale = flg ? MatDiagonalScale_Nest : 0;
961 
962  PetscFunctionReturn(0);
963 }
964 EXTERN_C_END
965 
966 #undef __FUNCT__
967 #define __FUNCT__ "MatNestSetVecType"
968 /*@C
969  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
970 
971  Not collective
972 
973  Input Parameters:
974 +  A  - nest matrix
975 -  vtype - type to use for creating vectors
976 
977  Notes:
978 
979  Level: developer
980 
981 .seealso: MatGetVecs()
982 @*/
983 PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
984 {
985   PetscErrorCode ierr;
986 
987   PetscFunctionBegin;
988   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
989   PetscFunctionReturn(0);
990 }
991 
992 EXTERN_C_BEGIN
993 #undef __FUNCT__
994 #define __FUNCT__ "MatNestSetSubMats_Nest"
995 PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
996 {
997   Mat_Nest       *s = (Mat_Nest*)A->data;
998   PetscInt       i,j,m,n,M,N;
999   PetscErrorCode ierr;
1000 
1001   PetscFunctionBegin;
1002   s->nr = nr;
1003   s->nc = nc;
1004 
1005   /* Create space for submatrices */
1006   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
1007   for (i=0; i<nr; i++) {
1008     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
1009   }
1010   for (i=0; i<nr; i++) {
1011     for (j=0; j<nc; j++) {
1012       s->m[i][j] = a[i*nc+j];
1013       if (a[i*nc+j]) {
1014         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1015       }
1016     }
1017   }
1018 
1019   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1020 
1021   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
1022   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
1023   for (i=0; i<nr; i++) s->row_len[i]=-1;
1024   for (j=0; j<nc; j++) s->col_len[j]=-1;
1025 
1026   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1027 
1028   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1029   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1030   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
1031   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1032   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1033   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
1034 
1035   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1036   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1037 
1038   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
1039   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
1040   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
1041   PetscFunctionReturn(0);
1042 }
1043 EXTERN_C_END
1044 
1045 #undef __FUNCT__
1046 #define __FUNCT__ "MatNestSetSubMats"
1047 /*@
1048    MatNestSetSubMats - Sets the nested submatrices
1049 
1050    Collective on Mat
1051 
1052    Input Parameter:
1053 +  N - nested matrix
1054 .  nr - number of nested row blocks
1055 .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1056 .  nc - number of nested column blocks
1057 .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1058 -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1059 
1060    Level: advanced
1061 
1062 .seealso: MatCreateNest(), MATNEST
1063 @*/
1064 PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1065 {
1066   PetscErrorCode ierr;
1067   PetscInt i;
1068 
1069   PetscFunctionBegin;
1070   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1071   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1072   if (nr && is_row) {
1073     PetscValidPointer(is_row,3);
1074     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1075   }
1076   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
1077   if (nc && is_col) {
1078     PetscValidPointer(is_col,5);
1079     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1080   }
1081   if (nr*nc) PetscValidPointer(a,6);
1082   ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr);
1083   PetscFunctionReturn(0);
1084 }
1085 
1086 /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1087 /*
1088   nprocessors = NP
1089   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
1090        proc 0: => (g_0,h_0,)
1091        proc 1: => (g_1,h_1,)
1092        ...
1093        proc nprocs-1: => (g_NP-1,h_NP-1,)
1094 
1095             proc 0:                      proc 1:                    proc nprocs-1:
1096     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
1097 
1098             proc 0:
1099     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
1100             proc 1:
1101     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
1102 
1103             proc NP-1:
1104     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
1105 */
1106 #undef __FUNCT__
1107 #define __FUNCT__ "MatSetUp_NestIS_Private"
1108 static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1109 {
1110   Mat_Nest       *vs = (Mat_Nest*)A->data;
1111   PetscInt       i,j,offset,n,nsum,bs;
1112   PetscErrorCode ierr;
1113   Mat            sub;
1114 
1115   PetscFunctionBegin;
1116   ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr);
1117   ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr);
1118   if (is_row) { /* valid IS is passed in */
1119     /* refs on is[] are incremeneted */
1120     for (i=0; i<vs->nr; i++) {
1121       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
1122       vs->isglobal.row[i] = is_row[i];
1123     }
1124   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
1125     nsum = 0;
1126     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
1127       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1128       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
1129       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
1130       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1131       nsum += n;
1132     }
1133     offset = 0;
1134     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1135     for (i=0; i<vs->nr; i++) {
1136       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1137       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
1138       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1139       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1140       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
1141       offset += n;
1142     }
1143   }
1144 
1145   if (is_col) { /* valid IS is passed in */
1146     /* refs on is[] are incremeneted */
1147     for (j=0; j<vs->nc; j++) {
1148       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
1149       vs->isglobal.col[j] = is_col[j];
1150     }
1151   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
1152     offset = A->cmap->rstart;
1153     nsum = 0;
1154     for (j=0; j<vs->nc; j++) {
1155       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1156       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
1157       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
1158       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1159       nsum += n;
1160     }
1161     offset = 0;
1162     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1163     for (j=0; j<vs->nc; j++) {
1164       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1165       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
1166       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1167       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1168       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
1169       offset += n;
1170     }
1171   }
1172 
1173   /* Set up the local ISs */
1174   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
1175   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
1176   for (i=0,offset=0; i<vs->nr; i++) {
1177     IS                     isloc;
1178     ISLocalToGlobalMapping rmap = PETSC_NULL;
1179     PetscInt               nlocal,bs;
1180     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1181     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);}
1182     if (rmap) {
1183       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1184       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1185       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1186       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1187     } else {
1188       nlocal = 0;
1189       isloc  = PETSC_NULL;
1190     }
1191     vs->islocal.row[i] = isloc;
1192     offset += nlocal;
1193   }
1194   for (i=0,offset=0; i<vs->nc; i++) {
1195     IS                     isloc;
1196     ISLocalToGlobalMapping cmap = PETSC_NULL;
1197     PetscInt               nlocal,bs;
1198     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1199     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);}
1200     if (cmap) {
1201       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1202       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1203       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1204       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1205     } else {
1206       nlocal = 0;
1207       isloc  = PETSC_NULL;
1208     }
1209     vs->islocal.col[i] = isloc;
1210     offset += nlocal;
1211   }
1212 
1213 #if defined(PETSC_USE_DEBUG)
1214   for (i=0; i<vs->nr; i++) {
1215     for (j=0; j<vs->nc; j++) {
1216       PetscInt m,n,M,N,mi,ni,Mi,Ni;
1217       Mat B = vs->m[i][j];
1218       if (!B) continue;
1219       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
1220       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
1221       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
1222       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
1223       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
1224       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1225       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);
1226       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);
1227     }
1228   }
1229 #endif
1230   PetscFunctionReturn(0);
1231 }
1232 
1233 #undef __FUNCT__
1234 #define __FUNCT__ "MatCreateNest"
1235 /*@
1236    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1237 
1238    Collective on Mat
1239 
1240    Input Parameter:
1241 +  comm - Communicator for the new Mat
1242 .  nr - number of nested row blocks
1243 .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1244 .  nc - number of nested column blocks
1245 .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1246 -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1247 
1248    Output Parameter:
1249 .  B - new matrix
1250 
1251    Level: advanced
1252 
1253 .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST
1254 @*/
1255 PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1256 {
1257   Mat            A;
1258   PetscErrorCode ierr;
1259 
1260   PetscFunctionBegin;
1261   *B = 0;
1262   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1263   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1264   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1265   *B = A;
1266   PetscFunctionReturn(0);
1267 }
1268 
1269 /*MC
1270   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1271 
1272   Level: intermediate
1273 
1274   Notes:
1275   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1276   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1277   It is usually used with DMComposite and DMGetMatrix()
1278 
1279 .seealso: MatCreate(), MatType, MatCreateNest()
1280 M*/
1281 EXTERN_C_BEGIN
1282 #undef __FUNCT__
1283 #define __FUNCT__ "MatCreate_Nest"
1284 PetscErrorCode MatCreate_Nest(Mat A)
1285 {
1286   Mat_Nest       *s;
1287   PetscErrorCode ierr;
1288 
1289   PetscFunctionBegin;
1290   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1291   A->data         = (void*)s;
1292   s->nr = s->nc   = -1;
1293   s->m            = PETSC_NULL;
1294 
1295   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1296   A->ops->mult                  = MatMult_Nest;
1297   A->ops->multadd               = 0;
1298   A->ops->multtranspose         = MatMultTranspose_Nest;
1299   A->ops->multtransposeadd      = 0;
1300   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1301   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1302   A->ops->zeroentries           = MatZeroEntries_Nest;
1303   A->ops->duplicate             = MatDuplicate_Nest;
1304   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1305   A->ops->destroy               = MatDestroy_Nest;
1306   A->ops->view                  = MatView_Nest;
1307   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1308   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1309   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1310   A->ops->getdiagonal           = MatGetDiagonal_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
1311   A->ops->diagonalscale         = MatDiagonalScale_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
1312 
1313   A->spptr        = 0;
1314   A->same_nonzero = PETSC_FALSE;
1315   A->assembled    = PETSC_FALSE;
1316 
1317   /* expose Nest api's */
1318   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1319   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C",   "MatNestSetSubMat_Nest",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1320   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1321   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1322   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1323   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1324 
1325   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1326   PetscFunctionReturn(0);
1327 }
1328 EXTERN_C_END
1329